fabi: Factor Analysis for Bicluster Acquisition: Laplace Prior...

Description Usage Arguments Details Value Author(s) References See Also Examples

View source: R/fabia.R

Description

fabi: R implementation of fabia, therefore it is slow.

Usage

1
fabi(X,p=13,alpha=0.01,cyc=500,spl=0,spz=0.5,center=2,norm=1,lap=1.0)

Arguments

X

the data matrix.

p

number of hidden factors = number of biclusters; default = 13.

alpha

sparseness loadings (0-1.0); default = 0.01.

cyc

number of iterations; default = 500.

spl

sparseness prior loadings (0 - 2.0); default = 0 (Laplace).

spz

sparseness factors (0.5-2.0); default = 0.5 (Laplace).

center

data centering: 1 (mean), 2 (median), > 2 (mode), 0 (no); default = 2.

norm

data normalization: 1 (0.75-0.25 quantile), >1 (var=1), 0 (no); default = 1.

lap

minimal value of the variational parameter; default = 1.0.

Details

Biclusters are found by sparse factor analysis where both the factors and the loadings are sparse.

Essentially the model is the sum of outer products of vectors:

X = ∑_{i=1}^{p} λ_i z_i^T + U

where the number of summands p is the number of biclusters. The matrix factorization is

X = L Z + U

Here λ_i are from R^n, z_i from R^l, L from R^{n \times p}, Z from R^{p \times l}, and X, U from R^{n \times l}.

If the nonzero components of the sparse vectors are grouped together then the outer product results in a matrix with a nonzero block and zeros elsewhere.

We recommend to normalize the components to variance one in order to make the signal and noise comparable across components.

The model selection is performed by a variational approach according to Girolami 2001 and Palmer et al. 2006.

We included a prior on the parameters and minimize a lower bound on the posterior of the parameters given the data. The update of the loadings includes an additive term which pushes the loadings toward zero (Gaussian prior leads to an multiplicative factor).

The code is implemented in R, therefore it is slow.

Value

object of the class Factorization. Containing LZ (estimated noise free data L Z), L (loadings L), Z (factors Z), U (noise X-LZ), center (centering vector), scaleData (scaling vector), X (centered and scaled data X), Psi (noise variance σ), lapla (variational parameter), avini (the information which the factor z_{ij} contains about x_j averaged over j) xavini (the information which the factor z_{j} contains about x_j) ini (for each j the information which the factor z_{ij} contains about x_j).

Author(s)

Sepp Hochreiter

References

S. Hochreiter et al., ‘FABIA: Factor Analysis for Bicluster Acquisition’, Bioinformatics 26(12):1520-1527, 2010. http://bioinformatics.oxfordjournals.org/cgi/content/abstract/btq227

Mark Girolami, ‘A Variational Method for Learning Sparse and Overcomplete Representations’, Neural Computation 13(11): 2517-2532, 2001.

J. Palmer, D. Wipf, K. Kreutz-Delgado, B. Rao, ‘Variational EM algorithms for non-Gaussian latent variable models’, Advances in Neural Information Processing Systems 18, pp. 1059-1066, 2006.

See Also

fabia, fabias, fabiap, spfabia, fabi, fabiasp, mfsc, nmfdiv, nmfeu, nmfsc, extractPlot, extractBic, plotBicluster, Factorization, projFuncPos, projFunc, estimateMode, makeFabiaData, makeFabiaDataBlocks, makeFabiaDataPos, makeFabiaDataBlocksPos, matrixImagePlot, fabiaDemo, fabiaVersion

Examples

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
#---------------
# TEST
#---------------

dat <- makeFabiaDataBlocks(n = 100,l= 50,p = 3,f1 = 5,f2 = 5,
  of1 = 5,of2 = 10,sd_noise = 3.0,sd_z_noise = 0.2,mean_z = 2.0,
  sd_z = 1.0,sd_l_noise = 0.2,mean_l = 3.0,sd_l = 1.0)

X <- dat[[1]]
Y <- dat[[2]]

resEx <- fabi(X,3,0.01,20)


## Not run: 
#---------------
# DEMO1
#---------------

dat <- makeFabiaDataBlocks(n = 1000,l= 100,p = 10,f1 = 5,f2 = 5,
  of1 = 5,of2 = 10,sd_noise = 3.0,sd_z_noise = 0.2,mean_z = 2.0,
  sd_z = 1.0,sd_l_noise = 0.2,mean_l = 3.0,sd_l = 1.0)

X <- dat[[1]]
Y <- dat[[2]]


resToy <- fabi(X,13,0.01,200)

extractPlot(resToy,ti="FABI",Y=Y)

#---------------
# DEMO2
#---------------

avail <- require(fabiaData)

if (!avail) {
    message("")
    message("")
    message("#####################################################")
    message("Package 'fabiaData' is not available: please install.")
    message("#####################################################")
} else {

data(Breast_A)

X <- as.matrix(XBreast)

resBreast <- fabi(X,5,0.1,200)

extractPlot(resBreast,ti="FABI Breast cancer(Veer)")

#sorting of predefined labels
CBreast
}

#---------------
# DEMO3
#---------------


avail <- require(fabiaData)

if (!avail) {
    message("")
    message("")
    message("#####################################################")
    message("Package 'fabiaData' is not available: please install.")
    message("#####################################################")
} else {


data(Multi_A)

X <- as.matrix(XMulti)

resMulti <- fabi(X,5,0.1,200)

extractPlot(resMulti,ti="FABI Multiple tissues(Su)")

#sorting of predefined labels
CMulti
}


#---------------
# DEMO4
#---------------


avail <- require(fabiaData)

if (!avail) {
    message("")
    message("")
    message("#####################################################")
    message("Package 'fabiaData' is not available: please install.")
    message("#####################################################")
} else {


data(DLBCL_B)

X <- as.matrix(XDLBCL)


resDLBCL <- fabi(X,5,0.1,200)

extractPlot(resDLBCL,ti="FABI Lymphoma(Rosenwald)")

#sorting of predefined labels
CDLBCL
}


## End(Not run)

fabia documentation built on Nov. 8, 2020, 8:09 p.m.