a dense or sparse explicit representation. x must be a
sparse explicit representation when a sparse kernel matrix should be
returned by the function (see parameter sparse).
y
a dense or sparse explicit representation. If x is dense,
y must be dense. If x is sparse, y must be sparse.
selx
a numeric or character vector for defining a subset of
x. Default=integer(0)
sely
a numeric or character vector for defining a subset of
y. Default=integer(0)
sparse
boolean indicating whether returned kernel matrix
should be sparse or dense. For value FALSE a dense kernel matrix
of class KernelMatrix is returned. If set to TRUE
the kernel matrix is returned as sparse matrix of class
dgCMatrix. In case of a symmetric matrix either the
lower triangular part or the full matrix can be returned. Please note that
a sparse kernel matrix currently can not be used for SVM based learning
in kebabs. Default=FALSE
triangular
boolean indicating whether just the lower triangular or
the full sparse matrix should be returned. This parameter is only relevant
for a sparse symmetric kernel matrix. Default=TRUE
diag
boolean indicating whether the diagonal should be included
in a sparse triangular matrix. This parameter is only relevant when
parameter sparse and triangular are set to TRUE.
Default=TRUE
lowerLimit
a numeric value for a similarity threshold. The
parameter is relevant for sparse kernel matrices only. If set to a value
larger than 0 only similarity values larger than this threshold will
be included in the sparse kernel matrix. Default=0
Value
linearKernel:
kernel matrix as class KernelMatrix or sparse
kernel matrix of class dgCMatrix
dependent on parameter sparse
J. Palme, S. Hochreiter, and U. Bodenhofer (2015) KeBABS: an R package
for kernel-based analysis of biological sequences.
Bioinformatics, 31(15):2574-2576, 2015.
DOI: 10.1093/bioinformatics/btv176.
Examples
## load sequence data and change sample names
data(TFBS)
names(enhancerFB) <- paste("S", 1:length(enhancerFB), sep="_")
## create the kernel object for dimers with normalization
speck <- spectrumKernel(k=5)
## generate sparse explicit representation
ers <- getExRep(enhancerFB, speck)
## compute dense kernel matrix (as currently used in SVM based learning)
km <- linearKernel(ers)
km[1:5, 1:5]
## compute sparse kernel matrix
## because it is symmetric just the lower diagonal
## is computed to save storage
km <- linearKernel(ers, sparse=TRUE)
km[1:5, 1:5]
## compute full sparse kernel matrix
km <- linearKernel(ers, sparse=TRUE, triangular=FALSE)
km[1:5, 1:5]
## compute triangular sparse kernel matrix without diagonal
km <- linearKernel(ers, sparse=TRUE, triangular=TRUE, diag=FALSE)
km[1:5, 1:5]
## plot histogram of similarity values
hist(as(km, "numeric"), breaks=30)
## compute sparse kernel matrix with similarities above 0.5 only
km <- linearKernel(ers, sparse=TRUE, lowerLimit=0.5)
km[1:5, 1:5]
Results
R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)
R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.
R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.
Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.
> library(kebabs)
Loading required package: Biostrings
Loading required package: BiocGenerics
Loading required package: parallel
Attaching package: 'BiocGenerics'
The following objects are masked from 'package:parallel':
clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
clusterExport, clusterMap, parApply, parCapply, parLapply,
parLapplyLB, parRapply, parSapply, parSapplyLB
The following objects are masked from 'package:stats':
IQR, mad, xtabs
The following objects are masked from 'package:base':
Filter, Find, Map, Position, Reduce, anyDuplicated, append,
as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
rbind, rownames, sapply, setdiff, sort, table, tapply, union,
unique, unsplit
Loading required package: S4Vectors
Loading required package: stats4
Attaching package: 'S4Vectors'
The following objects are masked from 'package:base':
colMeans, colSums, expand.grid, rowMeans, rowSums
Loading required package: IRanges
Loading required package: XVector
Loading required package: kernlab
Attaching package: 'kernlab'
The following object is masked from 'package:Biostrings':
type
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/kebabs/LinearKernel.Rd_%03d_medium.png", width=480, height=480)
> ### Name: linearKernel
> ### Title: Linear Kernel
> ### Aliases: linearKernel
> ### Keywords: kernel, linearKernel
>
> ### ** Examples
>
> ## load sequence data and change sample names
> data(TFBS)
> names(enhancerFB) <- paste("S", 1:length(enhancerFB), sep="_")
>
> ## create the kernel object for dimers with normalization
> speck <- spectrumKernel(k=5)
>
> ## generate sparse explicit representation
> ers <- getExRep(enhancerFB, speck)
>
> ## compute dense kernel matrix (as currently used in SVM based learning)
> km <- linearKernel(ers)
> km[1:5, 1:5]
An object of class "KernelMatrix"
S_1 S_2 S_3 S_4 S_5
S_1 1.0000000 0.4412064 0.4766687 0.5432537 0.5484077
S_2 0.4412064 1.0000000 0.4611439 0.5275230 0.4248695
S_3 0.4766687 0.4611439 1.0000000 0.5359756 0.4403322
S_4 0.5432537 0.5275230 0.5359756 1.0000000 0.5575523
S_5 0.5484077 0.4248695 0.4403322 0.5575523 1.0000000
>
> ## compute sparse kernel matrix
> ## because it is symmetric just the lower diagonal
> ## is computed to save storage
> km <- linearKernel(ers, sparse=TRUE)
> km[1:5, 1:5]
5 x 5 sparse Matrix of class "dgCMatrix"
S_1 S_2 S_3 S_4 S_5
S_1 1.0000000 . . . .
S_2 0.4412064 1.0000000 . . .
S_3 0.4766687 0.4611439 1.0000000 . .
S_4 0.5432537 0.5275230 0.5359756 1.0000000 .
S_5 0.5484077 0.4248695 0.4403322 0.5575523 1
>
> ## compute full sparse kernel matrix
> km <- linearKernel(ers, sparse=TRUE, triangular=FALSE)
> km[1:5, 1:5]
5 x 5 sparse Matrix of class "dgCMatrix"
S_1 S_2 S_3 S_4 S_5
S_1 1.0000000 0.4412064 0.4766687 0.5432537 0.5484077
S_2 0.4412064 1.0000000 0.4611439 0.5275230 0.4248695
S_3 0.4766687 0.4611439 1.0000000 0.5359756 0.4403322
S_4 0.5432537 0.5275230 0.5359756 1.0000000 0.5575523
S_5 0.5484077 0.4248695 0.4403322 0.5575523 1.0000000
>
> ## compute triangular sparse kernel matrix without diagonal
> km <- linearKernel(ers, sparse=TRUE, triangular=TRUE, diag=FALSE)
> km[1:5, 1:5]
5 x 5 sparse Matrix of class "dgCMatrix"
S_1 S_2 S_3 S_4 S_5
S_1 . . . . .
S_2 0.4412064 . . . .
S_3 0.4766687 0.4611439 . . .
S_4 0.5432537 0.5275230 0.5359756 . .
S_5 0.5484077 0.4248695 0.4403322 0.5575523 .
>
> ## plot histogram of similarity values
> hist(as(km, "numeric"), breaks=30)
>
> ## compute sparse kernel matrix with similarities above 0.5 only
> km <- linearKernel(ers, sparse=TRUE, lowerLimit=0.5)
> km[1:5, 1:5]
5 x 5 sparse Matrix of class "dgCMatrix"
S_1 S_2 S_3 S_4 S_5
S_1 1.0000000 . . . .
S_2 . 1.000000 . . .
S_3 . . 1.0000000 . .
S_4 0.5432537 0.527523 0.5359756 1.0000000 .
S_5 0.5484077 . . 0.5575523 1
>
>
>
>
>
> dev.off()
null device
1
>