Description Usage Arguments Details Value Table Frequency Constructor Author(s) Examples
Object that holds allele counts, genomic positions and map-bias for a set of SNPs
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 | alleleCounts(x, strand = "*", return.class = "list")
## S4 method for signature 'ASEset'
alleleCounts(x, strand = "*", return.class = "list")
alleleCounts(x, ...) <- value
## S4 replacement method for signature 'ASEset'
alleleCounts(x, strand = "*", return.class = "array", ...) <- value
mapBias(x, ...)
## S4 method for signature 'ASEset'
mapBias(x, return.class = "list")
fraction(x, ...)
## S4 method for signature 'ASEset'
fraction(
x,
strand = "*",
top.fraction.criteria = "maxcount",
verbose = FALSE,
...
)
arank(x, return.type = "names", return.class = "list", strand = "*", ...)
frequency(x, ...)
genotype(x, ...)
## S4 method for signature 'ASEset'
genotype(x, return.class = "matrix")
genotype(x) <- value
## S4 replacement method for signature 'ASEset'
genotype(x) <- value
countsPerSnp(x, ...)
## S4 method for signature 'ASEset'
countsPerSnp(x, return.class = "matrix", return.type = "mean", strand = "*")
countsPerSample(x, ...)
## S4 method for signature 'ASEset'
countsPerSample(x, return.class = "matrix", return.type = "mean", strand = "*")
phase(x, ...)
## S4 method for signature 'ASEset'
phase(x, return.class = "matrix")
phase(x) <- value
## S4 replacement method for signature 'ASEset'
phase(x) <- value
mapBias(x) <- value
## S4 replacement method for signature 'ASEset'
mapBias(x) <- value
refExist(x)
## S4 method for signature 'ASEset'
refExist(x)
ref(x)
## S4 method for signature 'ASEset'
ref(x)
ref(x) <- value
## S4 replacement method for signature 'ASEset,ANY'
ref(x) <- value
altExist(x)
## S4 method for signature 'ASEset'
altExist(x)
alt(x)
## S4 method for signature 'ASEset'
alt(x)
alt(x) <- value
## S4 replacement method for signature 'ASEset,ANY'
alt(x) <- value
aquals(x, ...)
## S4 method for signature 'ASEset'
aquals(x)
aquals(x) <- value
## S4 replacement method for signature 'ASEset'
aquals(x) <- value
maternalAllele(x, ...)
## S4 method for signature 'ASEset'
maternalAllele(x)
paternalAllele(x, ...)
## S4 method for signature 'ASEset'
paternalAllele(x)
|
x |
ASEset object |
strand |
which strand of '+', '-' or '*' |
return.class |
return 'list' or 'array' |
... |
additional arguments |
value |
replacement variable |
top.fraction.criteria |
'maxcount', 'ref' or 'phase' |
verbose |
makes function more talkative |
return.type |
return 'names', rank or 'counts' |
An ASEset object differs from a regular RangedSummarizedExperiment object in that the assays contains an array instead of matrix. This array has ranges on the rows, sampleNames on the columns and variants in the third dimension.
It is possible to use the commands barplot and locationplot on an ASEset
object see more details in barplot
and
locationplot
.
Three different alleleCount options are available. The simples one is the * option, and is for experiments where the strand information is not known e.g. non-stranded data. The unknown strand could also be for strand specific data when the aligner could not find any strand associated with the read, but this should normally not happen, and if it does probably having an extremely low mapping quality. Then there are an option too add plus and minus stranded data. When using this, it is essential to make sure that the RNA-seq experiment under analysis has in fact been created so that correct strand information was obtained. The most functions will by default have their strand argument set to '*'.
The phase information is stored by the convention of 'maternal chromosome|paternal chromosome', with 0 as reference allele and 1 as alternative allele. '|' when the phase is known and '/' when the phase is unknown. Internally the information will be stored as an three dimensional array, dim 1 for SNPs, dim 2 for Samples and dim 3 which is fixed and stores maternal chromosome, paternal chromosome and phased (1 equals TRUE).
An object of class ASEset containing location information and allele counts for a number of SNPs measured in a number of samples on various strand, as well as mapBias information. All data is stored in a manner similar to the SummarizedExperiment class.
table(...)
Arguments:
An ASEset object
that contains the
variants of interest
The generics for table does not easily allow more than one argument
so in respect to the different strand options, table
will
return a SimpleList with length 3, one element for each strand.
frequency(x, return.class = "list", strand = "*", threshold.count.sample = 15)
Arguments:
An ASEset object
that contains the
variants of interest
threshold.count.samples
if sample has fewer counts the function return NA.
ASEsetFromCountList(rowRanges, countListNonStranded = NULL, countListPlus = NULL, countListMinus = NULL, countListUnknown = NULL, colData = NULL, mapBiasExpMean = array(), verbose=FALSE, ...)
Arguments:
A GenomicRanges object
that contains the
variants of interest
A list
where each
entry is a matrix with allele counts as columns and sample counts as rows
A list
where each entry is a matrix with allele
counts as columns and sample counts as rows
A
list
where each entry is a matrix with allele counts as columns and
sample counts as rows
A list
where each
entry is a matrix with allele counts as columns and sample counts as rows
A DataFrame
object containing sample specific data
A 3D array
describing mapping bias. The SNPs
are in the 1st dimension, samples in the 2nd dimension and variants in the
3rd dimension.
Makes function more talkative
arguments passed on to SummarizedExperiment constructor
Jesper R. Gadin, Lasse Folkersen
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 | #make example countList
set.seed(42)
countListPlus <- list()
snps <- c('snp1','snp2','snp3','snp4','snp5')
for(snp in snps){
count<-matrix(rep(0,16),ncol=4,dimnames=list(
c('sample1','sample2','sample3','sample4'),
c('A','T','G','C')))
#insert random counts in two of the alleles
for(allele in sample(c('A','T','G','C'),2)){
count[,allele]<-as.integer(rnorm(4,mean=50,sd=10))
}
countListPlus[[snp]] <- count
}
#make example rowRanges
rowRanges <- GRanges(
seqnames = Rle(c('chr1', 'chr2', 'chr1', 'chr3', 'chr1')),
ranges = IRanges(1:5, width = 1, names = head(letters,5)),
snp = paste('snp',1:5,sep='')
)
#make example colData
colData <- DataFrame(Treatment=c('ChIP', 'Input','Input','ChIP'),
row.names=c('ind1','ind2','ind3','ind4'))
#make ASEset
a <- ASEsetFromCountList(rowRanges, countListPlus=countListPlus,
colData=colData)
#example phase matrix (simple form)
p1 <- matrix(sample(c(1,0),replace=TRUE, size=nrow(a)*ncol(a)),nrow=nrow(a), ncol(a))
p2 <- matrix(sample(c(1,0),replace=TRUE, size=nrow(a)*ncol(a)),nrow=nrow(a), ncol(a))
p <- matrix(paste(p1,sample(c("|","|","/"), size=nrow(a)*ncol(a), replace=TRUE), p2, sep=""),
nrow=nrow(a), ncol(a))
phase(a) <- p
#generate ASEset from array
snps <- 999
samples <-5
ar <-array(rep(unlist(lapply(1:snps,
function(x){(sample(c(TRUE,FALSE,TRUE,FALSE), size = 4))})), samples),
dim=c(4,snps,samples))
ar2 <- array(sample(50:300, 4*snps*samples,replace=TRUE), dim=c(4,snps,samples))
ar2[ar] <- 0
ar2 <- aperm(ar2, c(2, 3, 1))
dimnames(ar2) <- list(paste("snp",1:snps,sep=""),paste("sample",1:samples,sep=""),
c("A","C","G","T"))
gr <- GRanges(seqnames=c("chr2"), ranges=IRanges(start=1:dim(ar2)[1], width=1), strand="*")
a <- ASEsetFromArrays(gr, countsUnknown=ar2)
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.