AlignmentsExperimentSet-class | R Documentation |
Class AlignmentsExperimentSet
holds data from multiple alignments for
many experiments. Allows to examine alignments in great detail.
AlignmentsExperimentSet(...)
## S4 method for signature 'AlignmentsExperimentSet'
length(x)
## S4 method for signature 'AlignmentsExperimentSet'
fwdReads(x)
## S4 replacement method for signature 'AlignmentsExperimentSet'
fwdReads(x) <- value
## S4 method for signature 'AlignmentsExperimentSet'
rveReads(x)
## S4 replacement method for signature 'AlignmentsExperimentSet'
rveReads(x) <- value
## S4 method for signature 'AlignmentsExperimentSet'
fwdReadsType(x)
## S4 replacement method for signature 'AlignmentsExperimentSet'
fwdReadsType(x) <- value
## S4 method for signature 'AlignmentsExperimentSet'
rveReadsType(x)
## S4 replacement method for signature 'AlignmentsExperimentSet'
rveReadsType(x) <- value
## S4 method for signature 'AlignmentsExperimentSet'
unassignedData(x)
## S4 replacement method for signature 'AlignmentsExperimentSet'
unassignedData(x) <- value
## S4 method for signature 'AlignmentsExperimentSet'
readCounts(x)
## S4 replacement method for signature 'AlignmentsExperimentSet'
readCounts(x) <- value
## S4 method for signature 'AlignmentsExperimentSet'
experimentData(x)
## S4 replacement method for signature 'AlignmentsExperimentSet'
experimentData(x) <- value
## S4 method for signature 'AlignmentsExperimentSet'
barcodeData(x)
## S4 replacement method for signature 'AlignmentsExperimentSet'
barcodeData(x) <- value
## S4 method for signature 'AlignmentsExperimentSet'
unassignedCount(x)
## S4 method for signature 'AlignmentsExperimentSet'
assignedCount(x)
## S4 method for signature 'AlignmentsExperimentSet'
names(x)
## S4 method for signature 'AlignmentsExperimentSet'
c(x, ...)
## S4 method for signature 'AlignmentsExperimentSet,numeric,missing,missing'
x[i, j, ..., drop = TRUE]
## S3 method for class 'AlignmentsExperimentSet'
as.list(x, ...)
## S4 method for signature 'AlignmentsExperimentSet'
x$name
## S4 method for signature 'AlignmentsExperimentSet'
writeAlignments(x, file = "", aln_format = "txt")
## S4 method for signature 'AlignmentsExperimentSet'
lookupAlignment(x, ID, read_id = 1)
## S4 method for signature 'AlignmentsExperimentSet'
extractEvents(object, use_parallel = FALSE)
... |
pass any number of AlignmentsExperimentSet objects, make sure experiment IDs can be unique after merging |
x, object |
(AlignmentsExperimentSet) |
value |
Represents assignment values for setter methods. |
i, j, name, drop |
(numeric, missing, character, missing)
AlignmentsExperimentSet object can be subsetted using
names of the experiments eg. |
file |
(connection or string) Destination file. When empty, defaults to standard output. |
aln_format |
("txt" or "fasta") Specifies format of the file. |
ID |
(string) Experiment Identifier |
read_id |
(numeric) Read Identifier. Reads are sorted by frequency. Defaults to 1, most abundant read. |
use_parallel |
(boolean) When using |
depending on the function used
fwdReads,rveReads
(list) Named list where each element is of class
PairwiseAlignmentsSingleSubject
. Names correspond
to the experiment ID. Contains alignments of reads against amplicons.
fwdReadsType,rveReadsType
(list) Named list where each element is of logical vector, so far TRUE corresponds to HDR events. Names correspond to the experiment ID. Contains type of read - HDR/NHEJ.
readCounts
(list) Named list where each element is numeric vector that
describes how many reads are compressed into unique representation
before alignment in fwdReads
and/or rveReads
.
unassignedData
(data.frame) Contains reads that failed to be assigned to any of the experiments. Alignment of forward against reverse reads may give hint whether these reads are compromised in any way.
experimentData
(data.frame) Expands on configuration file and provides information about cut rates, frameshifts, PRIMER DIMER detection etc. Each row corresponds to experiment ID.
barcodeData
(data.frame) Information that is gathered on the barcode level is gathered in this data.frame, mainly quality filtering statistics.
Write out all alignments in "fasta" or "txt" format.:
writeAlignments(x, file = "", aln_format = "txt")
Write out human readable alignments for given experiment and read_id.:
lookupAlignment(x, ID, read_id = 1)
Coerce to data.frame
compatible with
GRanges
.:
as.data.frame(x)
exampleAlignments <- Biostrings::pairwiseAlignment(
Biostrings::DNAStringSet(c("ACTGACTG", "CGACGACG")), "ACGTACGTACGT")
new("AlignmentsExperimentSet",
fwdReads = list(ID_1 = exampleAlignments, ID_2 = exampleAlignments),
rveReads = list(ID_1 = exampleAlignments, ID_2 = exampleAlignments),
fwdReadsType = list(ID_1 = c(FALSE, FALSE), ID_2 = c(FALSE, FALSE)),
rveReadsType = list(ID_1 = c(FALSE, FALSE), ID_2 = c(FALSE, FALSE)),
readCounts = list(ID_1 = c(2, 20), ID_2 = c(30, 100)),
unassignedData = NULL,
experimentData = data.frame(ID = c("ID_1", "ID_2"),
Barcode = c("B1", "B1"),
whatever = c(50, 100)),
barcodeData = data.frame(Barcode = "B1", statistic1 = 100))
# Coercion
extractEvents(AlignmentsExperimentSet())
GenomicRanges::GRanges(extractEvents(AlignmentsExperimentSet()))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.