NanoStringRccSet-class: Class to Contain NanoString Expression Level Assays

NanoStringRccSet-classR Documentation

Class to Contain NanoString Expression Level Assays

Description

The NanoStringRccSet class extends the ExpressionSet class for NanoString Reporter Code Count (RCC) data.

Usage

NanoStringRccSet(assayData,
                 phenoData = annotatedDataFrameFrom(assayData, byrow = FALSE),
                 featureData = annotatedDataFrameFrom(assayData, byrow = TRUE),
                 experimentData = MIAME(),
                 annotation = character(),
                 protocolData = annotatedDataFrameFrom(assayData, byrow = FALSE),
                 dimLabels = c("GeneName", "SampleID"),
                 signatures = SignatureSet(),
                 design = NULL,
                 ...)

Arguments

assayData

A matrix or environment containing the RCCs.

phenoData

An AnnotatedDataFrame containing the phenotypic data.

featureData

An AnnotatedDataFrame containing columns "CodeClass", "GeneName", "Accession", "IsControl", and "ControlConc".

experimentData

An optional MIAME instance with meta-data about the experiment.

annotation

A character string for the "GeneRLF".

protocolData

An AnnotatedDataFrame containing columns "FileVersion", "SoftwareVersion", "SystemType", "SampleID", "SampleOwner", "SampleComments", "SampleDate", "SystemAPF", "AssayType", "LaneID", "FovCount", "FovCounted", "ScannerID", "StagePosition", "BindingDensity", "CartridgeID", and "CartridgeBarcode".

dimLabels

A character vector of length 2 that provides the column names to use as labels for the features and samples respectively in the autoplot method.

signatures

An optional SignatureSet object containing signature definitions.

design

An optional one-sided formula representing the experimental design based on columns from phenoData

...

Additional arguments for ExpressionSet.

Value

An S4 class containing NanoString Expression Level Assays

Accessing

In addition to the standard ExpressionSet accessor methods, NanoStringRccSet objects have the following:

sData(object)

extracts the data.frame containing the sample data, cbind(pData(object), pData(protocolData(object))).

svarLabels(object)

extracts the sample data column names, c(varLabels(object), varLabels(protocolData(object))).

dimLabels(object)

extracts the column names to use as labels for the features and samples in the autoplot method.

dimLabels(object) <- value

replaces the dimLabels of the object.

signatures(object)

extracts the SignatureSet of the object.

signatures(object) <- value

replaces the SignatureSet of the object.

signatureScores(object, elt = "exprs")

extracts the matrix of computed signature scores.

design(object)

extracts the one-sided formula representing the experimental design based on columns from phenoData.

design(object) <- value

replaces the one-sided formula representing the experimental design based on columns from phenoData.

setSignatureFuncs(object)

returns the signature functions.

setSignatureFuncs(object) <- value

replaces the signature functions.

setSignatureGroups(object) <- value

returns the signature groups.

setSignatureGroups(object) <- value

replaces the signature groups.

Summarizing

summary(object, MARGIN = 2L, GROUP = NULL, log2scale = TRUE, elt = "exprs", signatureScores = FALSE)

When signatureScores = FALSE, the marginal summaries of the elt assayData matrix along either the feature (MARGIN = 1) or sample (MARGIN = 2) dimension.

When signatureScores = TRUE, the marginal summaries of the elt signatureScores matrix along either the signature (MARGIN = 1) or sample (MARGIN = 2) dimension.

When log2scale = FALSE, the summary statistics are Mean, Standard Deviation, Skewness, Excess Kurtosis, Minimum, First Quartile, Median, Third Quartile, and Maximum.

When log2scale = TRUE, the summary statistics are Geometric Mean with thresholding at 0.5, Size Factor (2^(`MeanLog2` - mean(`MeanLog2`))), Mean of Log2 with thresholding at 0.5, Standard Deviation of Log2 with thresholding at 0.5, Minimum, First Quartile, Median, Third Quartile, and Maximum.

Subsetting

In addition to the standard ExpressionSet subsetting methods, NanoStringRccSet objects have the following:

subset(x, subset, select, ...)

Subset the feature and sample dimensions using the subset and select arguments respectively. The subset argument will be evaluated with respect to the featureData, while the select argument will be evaluated with respect to the phenoData and protocolData.

endogenousSubset(x, subset, select)

Extracts the endogenous barcode class feature subset of x with optional additional subsetting using subset and select.

housekeepingSubset(x, subset, select)

Extracts the housekeeping barcode class feature subset of x with optional additional subsetting using subset and select.

negativeControlSubset(x, subset, select)

Extracts the negative control barcode class feature subset of x with optional additional subsetting using subset and select.

positiveControlSubset(x, subset, select)

Extracts the positive control barcode class feature subset of x with optional additional subsetting using subset and select.

controlSubset(x, subset, select)

Extracts the feature subset representing the controls of x with optional additional subsetting using subset and select.

nonControlSubset(x, subset, select)

Extracts the feature subset representing the non-controls of x with optional additional subsetting using subset and select.

signatureSubset(x, subset, select)

Extracts the feature subset representing the genes in the signatures of x with optional additional subsetting using subset and select.

Looping

assayDataApply(X, MARGIN, FUN, ..., elt = "exprs")

Loop over the feature (MARGIN = 1) or sample (MARGIN = 2) dimension of assayDataElement(X, elt).

signatureScoresApply(X, MARGIN, FUN, ..., elt = "exprs")

Loop over the signature (MARGIN = 1) or sample (MARGIN = 2) dimension of signatureScores(X, elt).

esBy(X, GROUP, FUN, ..., simplify = TRUE)

Split X by GROUP column within featureData, phenoData, or protocolData and apply FUN to each partition.

Transforming

munge(data, mapping = update(design(data), exprs ~ .), extradata = NULL, elt = "exprs", ...)

munge argument data into a data.frame object for modeling and visualization using the mapping argument. Supplemental data can be specified using the extradata argument.

transform('_data', ...)

Similar to the transform generic in the base package, creates or modifies one or more assayData matrices based upon name = value pairs in .... The expressions in ... are appended to the preprocessing list in experimentData, which can be extracted using the preproc method.

Evaluating

with(data, expr, ...)

Evaluate expression expr with respect to assayData, featureData, phenoData, and protocolData; c(as.list(assayData(data)), fData(data), sData(data)).

Normalizing

normalize(object, type, fromElt = "exprs", toElt = "exprs_norm", ...)

Plotting

ggplot(data, mapping = aes(), ..., extradata = NULL, tooltip_digits = 4L, environment = parent.frame())

the NanoStringRccSet method for ggplot.

autoplot(object, type, log2scale = TRUE, elt = "exprs", index = 1L, geomParams = list(), tooltipDigits = 4L, heatmapGroup = NULL, ...)

Author(s)

Patrick Aboyoun

See Also

readNanoStringRccSet, writeNanoStringRccSet, ExpressionSet

Examples

# Create NanoStringRccSet from data files
datadir <- system.file("extdata", "3D_Bio_Example_Data",
                       package = "NanoStringNCTools")
rccs <- dir(datadir, pattern = "SKMEL.*\\.RCC$", full.names = TRUE)
rlf <- file.path(datadir, "3D_SolidTumor_Sig.rlf")
pheno <- file.path(datadir, "3D_SolidTumor_PhenoData.csv")
solidTumor <-
  readNanoStringRccSet(rccs, rlfFile = rlf, phenoDataFile = pheno)


# Create a deep copy of a NanoStringRccSet object
deepCopy <- NanoStringRccSet(solidTumor)
all.equal(solidTumor, deepCopy)
identical(solidTumor, deepCopy)


# Accessing sample data and column names
head(sData(solidTumor))
svarLabels(solidTumor)


# Set experimental design
design(solidTumor) <- ~ BRAFGenotype + Treatment
design(solidTumor)
munge(solidTumor)


# Marginal summarizing of NanoStringRccSet assayData matrices
head(summary(solidTumor, 1)) # Marginal summaries along features
head(summary(solidTumor, 2)) # Marginal summaries along samples


# Subsetting NanoStringRccSet objects
# Extract the positive controls for wildtype BRAF
dim(solidTumor)
dim(subset(solidTumor, CodeClass == "Positive", BRAFGenotype == "wt/wt"))

# Extract by barcode class
with(solidTumor, table(CodeClass))
with(endogenousSubset(solidTumor), table(CodeClass))
with(housekeepingSubset(solidTumor), table(CodeClass))
with(negativeControlSubset(solidTumor), table(CodeClass))
with(positiveControlSubset(solidTumor), table(CodeClass))
with(controlSubset(solidTumor), table(CodeClass))
with(nonControlSubset(solidTumor), table(CodeClass))


# Looping over NanoStringRccSet assayData matrices
log1pCoefVar <- function(x){
  x <- log1p(x)
  sd(x) / mean(x)
}

# Log1p Coefficient of Variation along Features
head(assayDataApply(solidTumor, 1, log1pCoefVar))

# Log1p Coefficient of Variation along Samples
head(assayDataApply(solidTumor, 2, log1pCoefVar))


# Transforming NanoSetRccSet assayData matrices
# Subtract max count from each sample
# Create log1p transformation of adjusted counts
thresh <- assayDataApply(negativeControlSubset(solidTumor), 2, max)
solidTumor2 <-
  transform(solidTumor,
            negCtrlZeroed = sweep(exprs, 2, thresh),
            log1p_negCtrlZeroed = log1p(pmax(negCtrlZeroed, 0)))
assayDataElementNames(solidTumor2)


# Evaluating expression using NanoStringRccSet data
meanLog1pExprs <-
  with(solidTumor,
       {
         means <- split(apply(exprs, 1, function(x) mean(log1p(x))), CodeClass)
         means <- means[order(sapply(means, median))]
         boxplot(means, horizontal = TRUE)
         means
       })

Nanostring-Biostats/NanoStringNCTools documentation built on April 19, 2024, 8:21 p.m.