estimateVBExpression: Estimate expression of transcripts using VB

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

View source: R/estimateVBExpression_call.R

Description

Estimates the expression of transcripts using Variational Bayes inference algorithm

Usage

1
2
3
estimateVBExpression (probFile, outFile, outputType=NULL, trInfoFile=NULL, 
      seed=NULL, samples=NULL, optLimit=1e-5, optMethod="FR", procN=4, 
      verbose=FALSE, veryVerbose=FALSE, pretend=FALSE)

Arguments

probFile

File with alignment probabilities produced by parseAlignment

outFile

Prefix for the output files.

outputType

Output type, possible values: theta, RPKM, counts. This is only relevant when the samples option is used. Default: theta.

trInfoFile

File containing transcript information. (Necessary for RPKM output)

seed

Sets the initial random seed for repeatable experiments.

samples

Number of samples to be generated from the posterior distribution. Default: no samples are generated.

verbose

Verbose output.

veryVerbose

Very verbose output.

procN

Maximum number of threads to be used. The program will not use more threads that there are MCMC chains.

Advanced options:

optLimit

The optimisation limit in terms of minimal gradient or change of bound.

optMethod

The optimisation method, use "FR", "HR", or "steepest".

pretend

Do not execute, only print out command line calls for the C++ version of the program.

Details

This function runs Variational Bayes algorithm to estimate the transcript expression. The input is the .prob file containing alignment probabilities which were produced by parseAlignment. Other optional input is the transcript information file specified by trInfoFile and again produced by parseAlignment.

It is much faster inference than MCMC which estimates mean expression equally well. However, the posterior is in form of Dirichlet distribution with underestimated variance. Use this method in cases when you are only interested in mean expression.

Value

.m_alphas

file containing mean relative expression of transcripts theta and parameters of the Dirichlet distribution. Please note the first line in the file corresponds to the noise transcript.

If option samples is used, the program also generates samples based of the outputType, the default would be file with extension ".VBtheta".

Author(s)

Peter Glaus

See Also

parseAlignment, estimateExpression

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
## Not run: 
setwd(system.file("extdata",package="BitSeq"))
parseAlignment( "data-c0b0.sam", outFile = "data-c0b0.prob", trSeqFile = "ensSelect1.fasta",
      trInfoFile = "data.tr", uniform = TRUE);

estimateVBExpression( probFile="data-c0b0.prob", outFile="data-c0b0-a", outputType="RPKM", 
      samples=1000, trInfoFile="data.tr", seed=47, verbose=TRUE)
estimateVBExpression( probFile="data-c0b0.prob", outFile="data-c0b0-b", trInfoFile="data.tr")
estimateVBExpression( probFile="data-c0b0.prob", outFile="data-c0b0-c", trInfoFile="data.tr", 
      optLimit=1e-6, optMethod = "HS", procN=12, veryVerbose=TRUE);

## End(Not run)

BitSeq documentation built on Nov. 8, 2020, 5:25 p.m.