Title: | NanoString GeoMx Tools |
---|---|
Description: | Tools for NanoString Technologies GeoMx Technology. Package provides functions for reading in DCC and PKC files based on an ExpressionSet derived object. Normalization and QC functions are also included. |
Authors: | Maddy Griswold [cre, aut], Nicole Ortogero [aut], Zhi Yang [aut], Ronalyn Vitancol [aut], David Henderson [aut] |
Maintainer: | Maddy Griswold <[email protected]> |
License: | MIT |
Version: | 3.11.0 |
Built: | 2024-10-30 07:29:30 UTC |
Source: | https://github.com/bioc/GeomxTools |
Aggregate probe counts to target level for feature data
aggregateCounts(object, FUN = ngeoMean)
aggregateCounts(object, FUN = ngeoMean)
object |
name of the NanoStringGeoMxSet object to aggregate |
FUN |
function to use for count aggregation |
a NanoStringGeoMxSet object with targets as features
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) targetGeoMxSet <- aggregateCounts(demoData[,1:10])
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) targetGeoMxSet <- aggregateCounts(demoData[,1:10])
Convert GeoMxSet Object to SeuratObject
## S3 method for class 'NanoStringGeoMxSet' as.Seurat( x, ident = NULL, normData = NULL, coordinates = NULL, forceRaw = FALSE, ... )
## S3 method for class 'NanoStringGeoMxSet' as.Seurat( x, ident = NULL, normData = NULL, coordinates = NULL, forceRaw = FALSE, ... )
x |
An object to convert to class |
ident |
column in GeoMxSet segmentProperties to set as Seurat object's identity class |
normData |
assay containing normalized data |
coordinates |
X and Y coordinates of each ROI, format: c(X,Y) |
forceRaw |
should raw data be forced into SeuratObject, not recommended |
... |
Arguments passed to other methods |
SeuratObject containing GeoMx data
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) target_demoData <- aggregateCounts(demoData[1:1000,1:10]) target_demoData <- normalize(target_demoData, "quant") seurat_demoData <- as.Seurat(target_demoData, ident = "cell_line", normData = "exprs_norm", forceRaw = FALSE)
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) target_demoData <- aggregateCounts(demoData[1:1000,1:10]) target_demoData <- normalize(target_demoData, "quant") seurat_demoData <- as.Seurat(target_demoData, ident = "cell_line", normData = "exprs_norm", forceRaw = FALSE)
Convert Object to SpatialExperiment
Convert GeoMxSet Object to SpatialExperiment
as.SpatialExperiment(x, ...) ## S3 method for class 'NanoStringGeoMxSet' as.SpatialExperiment( x, normData = NULL, coordinates = NULL, forceRaw = FALSE, ... )
as.SpatialExperiment(x, ...) ## S3 method for class 'NanoStringGeoMxSet' as.SpatialExperiment( x, normData = NULL, coordinates = NULL, forceRaw = FALSE, ... )
x |
GeoMxSet object to convert |
... |
Arguments passed to other methods |
normData |
assay containing normalized data |
coordinates |
X and Y coordinates of each ROI, format: c(X,Y) |
forceRaw |
should raw data be forced into SpatialExperiment, not recommended |
SpatialExperiment containing GeoMx data
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) target_demoData <- aggregateCounts(demoData[1:1000,1:10]) target_demoData <- normalize(target_demoData, "quant") seurat_demoData <- as.SpatialExperiment(target_demoData, normData = "exprs_norm", forceRaw = FALSE)
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) target_demoData <- aggregateCounts(demoData[1:1000,1:10]) target_demoData <- normalize(target_demoData, "quant") seurat_demoData <- as.SpatialExperiment(target_demoData, normData = "exprs_norm", forceRaw = FALSE)
Check QC Flags in the GeoMxSet and removes the probe or sample from the object
checkQCFlags(object, ...)
checkQCFlags(object, ...)
object |
name of the NanoStringGeoMxSet object to check the QC Flags |
... |
for other arguments |
a NanoStringGeoMxSet object probes and samples failing QC removed
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools" ) demoData <- readRDS(file.path(datadir, "/demoData.rds")) QCobject <- checkQCFlags(demoData)
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools" ) demoData <- readRDS(file.path(datadir, "/demoData.rds")) QCobject <- checkQCFlags(demoData)
checkQCFlags
## S4 method for signature 'NanoStringGeoMxSet' checkQCFlags(object, removeLowLocalOutliers = FALSE, ...)
## S4 method for signature 'NanoStringGeoMxSet' checkQCFlags(object, removeLowLocalOutliers = FALSE, ...)
object |
name of the NanoStringGeoMxSet object to check the QC Flags |
removeLowLocalOutliers |
logical, if TRUE it sets outlier counts to zero, default is FALSE, |
... |
optional arguments |
NanoStringGeoMxSet
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools" ) demoData <- readRDS(file.path(datadir, "/demoData.rds")) QCobject <- checkQCFlags(demoData)
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools" ) demoData <- readRDS(file.path(datadir, "/demoData.rds")) QCobject <- checkQCFlags(demoData)
Check if extra PKCs are given based on probes in config file
compareToConfig(config, pkcProbes, pkcHeader)
compareToConfig(config, pkcProbes, pkcHeader)
config |
file path to config file |
pkcProbes |
probe information from readPKCFile |
pkcHeader |
pkc metadata from readPKCFile |
For use with protein data ONLY.
Generate normalization factors for protein data to determine the best normalization method
computeNormalizationFactors( object, igg.names = NULL, hk.names = NULL, area = NULL, nuclei = NULL )
computeNormalizationFactors( object, igg.names = NULL, hk.names = NULL, area = NULL, nuclei = NULL )
object |
name of the object class to subset
|
igg.names |
names of IgGs, if NULL IgGs will be detected automatically |
hk.names |
names of HK, if NULL HK will be detected automatically |
area |
name of area column in annotation sheet, optional |
nuclei |
name of nuclei column in annotation sheet, optional |
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) normfactors <- computeNormalizationFactors(object = proteinData) normfactors_withAreaNuclei <- computeNormalizationFactors(object = proteinData, area = "AOI.Size.um2", nuclei = "Nuclei.Counts")
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) normfactors <- computeNormalizationFactors(object = proteinData) normfactors_withAreaNuclei <- computeNormalizationFactors(object = proteinData, area = "AOI.Size.um2", nuclei = "Nuclei.Counts")
assDataElement
was shifted by oneAccessor to check if "exprs" assDataElement
was shifted by one
countsShiftedByOne(object)
countsShiftedByOne(object)
object |
name of the NanoStringGeoMxSet object |
boolean indicating if counts in default matrix were shifted by one
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) countsShiftedByOne(demoData)
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) countsShiftedByOne(demoData)
Return the House Keeper positive controls for protein
hkNames(object)
hkNames(object)
object |
name of the NanoStringGeoMxSet object |
names of HKs
Return the IgG negative controls for protein
iggNames(object)
iggNames(object)
object |
name of the NanoStringGeoMxSet object |
names of IgGs
Get take the log of a numeric vector
logtBase(x, thresh = 0.5, base = 2)
logtBase(x, thresh = 0.5, base = 2)
x |
numeric vector |
thresh |
minimum numeric value greater than 0 to have in vector |
base |
numeric value indicating base to log with |
numeric vector with logged values
logtBase(c(0, 1, 2, 2), thresh=0.1, base=10)
logtBase(c(0, 1, 2, 2), thresh=0.1, base=10)
Run a mixed model on GeoMxSet
mixedModelDE( object, elt = "exprs", modelFormula = NULL, groupVar = "group", nCores = 1, multiCore = TRUE, pAdjust = "BY", pairwise = TRUE )
mixedModelDE( object, elt = "exprs", modelFormula = NULL, groupVar = "group", nCores = 1, multiCore = TRUE, pAdjust = "BY", pairwise = TRUE )
object |
name of the object class to perform QC on
|
elt |
assayDataElement of the geoMxSet object to run the DE on |
modelFormula |
formula used in DE, if null, the design(object) is used |
groupVar |
= "group", sample annotation to group the data for comparing means |
nCores |
= 1, number of cores to use, set to 1 if running in serial mode |
multiCore |
= TRUE, set to TRUE to use multiCore, FALSE to run in cluster mode |
pAdjust |
= "BY" method for p-value adjustment |
pairwise |
boolean to calculate least-square means pairwise differences |
mixed model output list
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) target_demoData <- aggregateCounts(demoData) target_demoData <- normalize(target_demoData, norm_method="quant") pData(target_demoData)[["slide"]] <- factor(pData(target_demoData)[["slide name"]]) protocolData(target_demoData)[["pool_rep"]] <- factor(protocolData(target_demoData)[["pool_rep"]]) mixedOutmc <- mixedModelDE(target_demoData, elt = "exprs_norm", modelFormula = ~ pool_rep + (1 | slide), groupVar = "pool_rep", nCores = 2, multiCore = TRUE, pAdjust = NULL )
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) target_demoData <- aggregateCounts(demoData) target_demoData <- normalize(target_demoData, norm_method="quant") pData(target_demoData)[["slide"]] <- factor(pData(target_demoData)[["slide name"]]) protocolData(target_demoData)[["pool_rep"]] <- factor(protocolData(target_demoData)[["pool_rep"]]) mixedOutmc <- mixedModelDE(target_demoData, elt = "exprs_norm", modelFormula = ~ pool_rep + (1 | slide), groupVar = "pool_rep", nCores = 2, multiCore = TRUE, pAdjust = NULL )
The NanoStringGeoMxSet
class extends the
ExpressionSet
class for NanoString GeoMx Digital Count Conversion
(DCC) data.
NanoStringGeoMxSet(assayData, phenoData=Biobase::annotatedDataFrameFrom(assayData, byrow=FALSE), featureData=Biobase::annotatedDataFrameFrom(assayData, byrow=TRUE), experimentData=Biobase::MIAME(), annotation=character(), protocolData=Biobase::annotatedDataFrameFrom(assayData, byrow=FALSE), dimLabels=c("TargetName", "SampleID"), signatures=SignatureSet(), design=NULL, featureType="Probe", analyte="RNA", ...)
NanoStringGeoMxSet(assayData, phenoData=Biobase::annotatedDataFrameFrom(assayData, byrow=FALSE), featureData=Biobase::annotatedDataFrameFrom(assayData, byrow=TRUE), experimentData=Biobase::MIAME(), annotation=character(), protocolData=Biobase::annotatedDataFrameFrom(assayData, byrow=FALSE), dimLabels=c("TargetName", "SampleID"), signatures=SignatureSet(), design=NULL, featureType="Probe", analyte="RNA", ...)
assayData |
A |
phenoData |
An |
featureData |
An |
experimentData |
An optional |
annotation |
A |
protocolData |
An |
dimLabels |
A |
signatures |
An optional |
design |
An optional one-sided formula representing the experimental
design based on columns from |
featureType |
A |
analyte |
A |
... |
Additional arguments for |
An S4 class containing data from a NanoString GeoMx experiment
Objects can be updated to current version using updateGeoMxSet(object)
In addition to the standard ExpressionSet
accessor
methods, NanoStringGeoMxSet
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.
dimLabels(object) <- value
replaces the dimLabels
of
the object
.
featureType(object)
extracts the featureType
of the object
.
featureType(object) <- value
replaces the featureType
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
.
signatureGroups(object)
extract the groups of SignatureSet
.
signatureGroups(object) <- value
replaces the groups of SignatureSet
.
analyte(object)
extracts the analyte
of the object
.
Zhi Yang & Nicole Ortogero
readNanoStringGeoMxSet
,
ExpressionSet
# Create NanoStringGeoMxSet from data files datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") dccFiles <- dir(datadir, pattern=".dcc$", full.names=TRUE) pkc <- unzip(zipfile = file.path(datadir, "/pkcs.zip")) sampleAnnotationFile <- file.path(datadir, "annotations.xlsx") dccFileColumn <- "Sample_ID" dccSet <- readNanoStringGeoMxSet(dccFiles=dccFiles, pkcFiles=pkc, phenoDataFile=sampleAnnotationFile, phenoDataSheet="CW005", phenoDataDccColName=dccFileColumn, protocolDataColNames=c("aoi", "cell_line", "roi_rep", "pool_rep", "slide_rep"), experimentDataColNames="panel", phenoDataColPrefix="") # Accessing sample data and column names head(sData(dccSet)) svarLabels(dccSet) featureType(dccSet) analyte(dccSet) # Accessing number of samples and features dim(dccSet)
# Create NanoStringGeoMxSet from data files datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") dccFiles <- dir(datadir, pattern=".dcc$", full.names=TRUE) pkc <- unzip(zipfile = file.path(datadir, "/pkcs.zip")) sampleAnnotationFile <- file.path(datadir, "annotations.xlsx") dccFileColumn <- "Sample_ID" dccSet <- readNanoStringGeoMxSet(dccFiles=dccFiles, pkcFiles=pkc, phenoDataFile=sampleAnnotationFile, phenoDataSheet="CW005", phenoDataDccColName=dccFileColumn, protocolDataColNames=c("aoi", "cell_line", "roi_rep", "pool_rep", "slide_rep"), experimentDataColNames="panel", phenoDataColPrefix="") # Accessing sample data and column names head(sData(dccSet)) svarLabels(dccSet) featureType(dccSet) analyte(dccSet) # Accessing number of samples and features dim(dccSet)
Get the geometric mean of a vector
ngeoMean(x, thresh = 0.5)
ngeoMean(x, thresh = 0.5)
x |
numeric vector |
thresh |
minimum numeric value greater than 0 to have in vector |
numeric geometric mean of vector
ngeoMean(c(0, 1, 2, 2), thresh=0.1)
ngeoMean(c(0, 1, 2, 2), thresh=0.1)
Get the geometric standard deviation of a vector
ngeoSD(x, thresh = 0.5)
ngeoSD(x, thresh = 0.5)
x |
numeric vector |
thresh |
minimum numeric value greater than 0 to have in vector |
numeric geometric standard deviation of vector
ngeoSD(c(0, 1, 2, 2), thresh=0.1)
ngeoSD(c(0, 1, 2, 2), thresh=0.1)
normalize GeoMxSet using different normalization methods
## S4 method for signature 'NanoStringGeoMxSet' normalize( object, norm_method = c("quant", "neg", "hk", "subtractBackground"), fromElt = "exprs", toElt = "exprs_norm", housekeepers = HOUSEKEEPERS, ... )
## S4 method for signature 'NanoStringGeoMxSet' normalize( object, norm_method = c("quant", "neg", "hk", "subtractBackground"), fromElt = "exprs", toElt = "exprs_norm", housekeepers = HOUSEKEEPERS, ... )
object |
name of the object class to perform normalization on |
norm_method |
the normalization method to be applied on the object |
fromElt |
name of the assayDataElement to normalize |
toElt |
name of the assayDataElement to store normalized values |
housekeepers |
optional vector of housekeeper target names |
... |
optional arguments |
a NanoStringGeoMxSet object with normalized counts and normalized factors
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools" ) demoData <- readRDS(file.path(datadir, "/demoData.rds")) norm_object <- normalize(demoData[1:1000,1:10])
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package = "GeomxTools" ) demoData <- readRDS(file.path(datadir, "/demoData.rds")) norm_object <- normalize(demoData[1:1000,1:10])
Upper panels are the concordance plot. Lower panels are the standard deviation of the log2-ratios between the targets
plotConcordance(targetList, object, plotFactor)
plotConcordance(targetList, object, plotFactor)
targetList |
names of targets to plot concordance, normally IgGs. |
object |
name of the object class to subset
|
plotFactor |
segment factor to color the plot by |
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) igg.names <- iggNames(proteinData) protSegTypeFig <- plotConcordance(targetList = igg.names, object = proteinData, plotFactor = "Segment_Type") protSegTypeFig
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) igg.names <- iggNames(proteinData) protSegTypeFig <- plotConcordance(targetList = igg.names, object = proteinData, plotFactor = "Segment_Type") protSegTypeFig
For use with protein data ONLY.
Upper panels are the concordance plot. Lower panels are the standard deviation of the log2-ratios between the normalization factors
plotNormFactorConcordance(object, plotFactor, normfactors = NULL)
plotNormFactorConcordance(object, plotFactor, normfactors = NULL)
object |
name of the object class to subset
|
plotFactor |
segment factor to color the plot by |
normfactors |
normalization factors from computeNormalizationFactors(). If NULL these are calculated automatically. |
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) normConcord <- plotNormFactorConcordance(object = proteinData, plotFactor = "Segment_Type") normConcord
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) normConcord <- plotNormFactorConcordance(object = proteinData, plotFactor = "Segment_Type") normConcord
For use with protein data ONLY.
qcProteinSignal(object, neg.names = NULL)
qcProteinSignal(object, neg.names = NULL)
object |
name of the object class to subset
|
neg.names |
names of IgGs, if NULL IgGs will be detected automatically |
figure function
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) igg.names <- iggNames(proteinData) qcFig <- qcProteinSignal(object = proteinData, neg.names = igg.names) qcFig()
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) igg.names <- iggNames(proteinData) qcFig <- qcProteinSignal(object = proteinData, neg.names = igg.names) qcFig()
For use with protein data ONLY.
qcProteinSignalNames(object, neg.names)
qcProteinSignalNames(object, neg.names)
object |
name of the object class to subset
|
neg.names |
names of IgGs, if NULL IgGs will be detected automatically |
protein names in increasing SNR order
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) igg.names <- iggNames(proteinData) proteinOrder <- qcProteinSignalNames(object = proteinData, neg.names = igg.names)
proteinData <- readRDS(file= system.file("extdata","DSP_Proteogenomics_Example_Data", "proteinData.rds", package = "GeomxTools")) igg.names <- iggNames(proteinData) proteinOrder <- qcProteinSignalNames(object = proteinData, neg.names = igg.names)
Read a NanoString GeoMx Digital Count Conversion (DCC) file.
readDccFile(file)
readDccFile(file)
file |
A character string containing the path to the DCC file. |
A list object with two elements:
"Header" |
a |
"Code_Summary" |
a |
Zhi Yang & Nicole Ortogero
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") dccFiles <- dir(datadir, pattern=".dcc$", full.names=TRUE) dccData <- sapply(dccFiles[1:10], readDccFile, simplify = FALSE)
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") dccFiles <- dir(datadir, pattern=".dcc$", full.names=TRUE) dccData <- sapply(dccFiles[1:10], readDccFile, simplify = FALSE)
Create an instance of class NanoStringGeoMxSet
by reading
data from NanoString GeoMx Digital Count Conversion (DCC) data.
readNanoStringGeoMxSet(dccFiles, pkcFiles, phenoDataFile, phenoDataSheet, phenoDataDccColName = "Sample_ID", phenoDataColPrefix = "", protocolDataColNames = NULL, experimentDataColNames = NULL, configFile = NULL, analyte = "RNA", defaultPKCVersions = NULL, ...)
readNanoStringGeoMxSet(dccFiles, pkcFiles, phenoDataFile, phenoDataSheet, phenoDataDccColName = "Sample_ID", phenoDataColPrefix = "", protocolDataColNames = NULL, experimentDataColNames = NULL, configFile = NULL, analyte = "RNA", defaultPKCVersions = NULL, ...)
dccFiles |
A character vector containing the paths to the DCC files. |
pkcFiles |
A character vector representing the path to the corresponding PKC file. |
phenoDataFile |
Character string representing the path to the corresponding phenotypic excel data file. It is recommended to use the Lab Worksheet in the exact order samples are provided in. |
phenoDataSheet |
Character string representing the excel sheet name containing the phenotypic data. |
phenoDataDccColName |
Character string identifying unique sample identifier
column in |
phenoDataColPrefix |
An optional prefix to add to the phenoData column names to distinguish them from the names of assayData matrices, featureData columns, and protocolData columns. |
protocolDataColNames |
Character list of column names from |
experimentDataColNames |
Character list of column names from |
configFile |
An optional character string representing the path to the corresponding config file. This is used to ensure the only the correct PKC files are added |
analyte |
GeoMxSet objects can only hold one analyte at a time. For studies with multiple analytes, which one should be read in? Options: RNA (default) and Protein |
defaultPKCVersions |
Optional list of pkc file names to use as default if more than one pkc version of each module is provided. |
... |
Optional parameters to pass to |
An instance of the NanoStringGeoMxSet
class.
Zhi Yang & Nicole Ortogero
# Data file paths datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") dccFiles <- dir(datadir, pattern=".dcc$", full.names=TRUE) pkc <- unzip(zipfile = file.path(datadir, "/pkcs.zip")) sampleAnnotationFile <- file.path(datadir, "annotations.xlsx") dccFileColumn <- "Sample_ID" dccSet <- readNanoStringGeoMxSet(dccFiles=dccFiles[1:10], pkcFiles=pkc, phenoDataFile=sampleAnnotationFile, phenoDataSheet="CW005", phenoDataDccColName=dccFileColumn, protocolDataColNames=c("aoi", "cell_line", "roi_rep", "pool_rep", "slide_rep"), experimentDataColNames="panel", phenoDataColPrefix="") # All data dccSet <- readNanoStringGeoMxSet(dccFiles, pkcFile = pkc, phenoDataFile = sampleAnnotationFile, phenoDataSheet="CW005") varLabels(dccSet) # All data with phenoData prefix dccSetPhenoPrefix <- readNanoStringGeoMxSet(dccFiles, pkcFile = pkc, phenoDataFile = sampleAnnotationFile, phenoDataSheet="CW005", phenoDataColPrefix = "PHENO_") varLabels(dccSetPhenoPrefix)
# Data file paths datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") dccFiles <- dir(datadir, pattern=".dcc$", full.names=TRUE) pkc <- unzip(zipfile = file.path(datadir, "/pkcs.zip")) sampleAnnotationFile <- file.path(datadir, "annotations.xlsx") dccFileColumn <- "Sample_ID" dccSet <- readNanoStringGeoMxSet(dccFiles=dccFiles[1:10], pkcFiles=pkc, phenoDataFile=sampleAnnotationFile, phenoDataSheet="CW005", phenoDataDccColName=dccFileColumn, protocolDataColNames=c("aoi", "cell_line", "roi_rep", "pool_rep", "slide_rep"), experimentDataColNames="panel", phenoDataColPrefix="") # All data dccSet <- readNanoStringGeoMxSet(dccFiles, pkcFile = pkc, phenoDataFile = sampleAnnotationFile, phenoDataSheet="CW005") varLabels(dccSet) # All data with phenoData prefix dccSetPhenoPrefix <- readNanoStringGeoMxSet(dccFiles, pkcFile = pkc, phenoDataFile = sampleAnnotationFile, phenoDataSheet="CW005", phenoDataColPrefix = "PHENO_") varLabels(dccSetPhenoPrefix)
Read a NanoString Probe Kit Configuration (PKC) file.
readPKCFile(file, default_pkc_vers=NULL)
readPKCFile(file, default_pkc_vers=NULL)
file |
A character string containing the path to the PKC file. |
default_pkc_vers |
Optional list of pkc file names to use as default if more than one pkc version of each module is provided. |
An instance of the DataFrame
class containing columns:
"RTS_ID" |
unique probe ID |
"TargetName" |
target or gene name |
"Module" |
PKC name |
"Negative" |
negative probe |
... |
additional columns |
Zhi Yang & Nicole Ortogero
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") pkc <- unzip(zipfile = file.path(datadir, "/pkcs.zip")) PKCData <- readPKCFile(pkc)
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") pkc <- unzip(zipfile = file.path(datadir, "/pkcs.zip")) PKCData <- readPKCFile(pkc)
Add background QC flags to NanoStringGeoMxSet object protocol data
setBackgroundQCFlags(object, qcCutoffs = DEFAULTS)
setBackgroundQCFlags(object, qcCutoffs = DEFAULTS)
object |
name of the NanoStringGeoMxSet object to perform QC on |
qcCutoffs |
a list of qc cutoffs to use
|
NanoStringGeoMxSet
object with QCFlags
data frame
appended to protocolData
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setBackgroundQCFlags(demoData[,1:10], qcCutoffs=list(minNegativeCount=10, maxNTCCount=60))
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setBackgroundQCFlags(demoData[,1:10], qcCutoffs=list(minNegativeCount=10, maxNTCCount=60))
Add probe QC flags to NanoStringGeoMxSet object feature data
setBioProbeQCFlags(object, qcCutoffs = DEFAULTS, removeLocalOutliers = TRUE)
setBioProbeQCFlags(object, qcCutoffs = DEFAULTS, removeLocalOutliers = TRUE)
object |
name of the NanoStringGeoMxSet object to perform QC on |
qcCutoffs |
a list of qc cutoffs to use
|
removeLocalOutliers |
boolean to determine if
local outliers should be excluded from |
NanoStringGeoMxSet
object with QCFlags
data frame
appended to protocolData
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) demoData <- shiftCountsOne(demoData, elt="exprs", useDALogic=TRUE) setBioProbeQCFlags(demoData[,1:10], qcCutoffs=list(minProbeRatio=0.1, percentFailGrubbs=20), removeLocalOutliers=TRUE)
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) demoData <- shiftCountsOne(demoData, elt="exprs", useDALogic=TRUE) setBioProbeQCFlags(demoData[,1:10], qcCutoffs=list(minProbeRatio=0.1, percentFailGrubbs=20), removeLocalOutliers=TRUE)
Add GeoMx segment QC flags to NanoStringGeoMxSet object protocol data
setGeoMxQCFlags(object, qcCutoffs = DEFAULTS)
setGeoMxQCFlags(object, qcCutoffs = DEFAULTS)
object |
name of the NanoStringGeoMxSet object to perform QC on |
qcCutoffs |
a list of qc cutoffs to use
|
NanoStringGeoMxSet
object with QCFlags
data frame
appended to protocolData
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setGeoMxQCFlags(demoData[,1:10], qcCutoffs=list(minNuclei=16000, minArea=20))
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setGeoMxQCFlags(demoData[,1:10], qcCutoffs=list(minNuclei=16000, minArea=20))
Add QC flags to feature and protocol data simultaneously
## S4 method for signature 'NanoStringGeoMxSet' setQCFlags(object, qcCutoffs = DEFAULTS, ...)
## S4 method for signature 'NanoStringGeoMxSet' setQCFlags(object, qcCutoffs = DEFAULTS, ...)
object |
name of the object class to perform QC on
|
qcCutoffs |
list of cutoffs and thresholds to use for QC |
... |
optional parameters to pass |
the object that QC was performed on
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setQCFlags(demoData[,1:10])
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setQCFlags(demoData[,1:10])
Add segment QC flags to protocol data
setSegmentQCFlags(object, qcCutoffs = DEFAULTS)
setSegmentQCFlags(object, qcCutoffs = DEFAULTS)
object |
name of the object class to perform QC on
|
qcCutoffs |
list of cutoffs and thresholds to use for QC |
NanoStringGeoMxSet
object with QCFlags
data frame
appended to protocolData
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setSegmentQCFlags(demoData[,1:10], qcCutoffs=list(minSegmentReads=1000, percentAligned=80, percentSaturation=50, minNegativeCount=10, maxNTCCount=60, minNuclei=16000, minArea=20))
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setSegmentQCFlags(demoData[,1:10], qcCutoffs=list(minSegmentReads=1000, percentAligned=80, percentSaturation=50, minNegativeCount=10, maxNTCCount=60, minNuclei=16000, minArea=20))
Add sequencing QC flags to NanoStringGeoMxSet object protocol data
setSeqQCFlags(object, qcCutoffs = DEFAULTS)
setSeqQCFlags(object, qcCutoffs = DEFAULTS)
object |
name of the NanoStringGeoMxSet object to perform QC on |
qcCutoffs |
a list of qc cutoffs to use
|
NanoStringGeoMxSet
object with QCFlags
data frame
appended to protocolData
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setSeqQCFlags(demoData[,1:10], qcCutoffs=list(minSegmentReads=1000, percentAligned=80, percentSaturation=50))
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) setSeqQCFlags(demoData[,1:10], qcCutoffs=list(minSegmentReads=1000, percentAligned=80, percentSaturation=50))
Add one to all counts in an expression matrix
shiftCountsOne(object, elt = "exprs", useDALogic = FALSE)
shiftCountsOne(object, elt = "exprs", useDALogic = FALSE)
object |
name of the NanoStringGeoMxSet object |
elt |
expression matrix element in |
useDALogic |
boolean to use the same logic in DA (impute 0s to 1s) or set to FALSE to shift all counts by 1 |
object of NanoStringGeoMxSet class
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) shiftCountsOne(demoData)
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) shiftCountsOne(demoData)
Calculate negative probe summary stats
summarizeNegatives(object, functionList = c())
summarizeNegatives(object, functionList = c())
object |
name of the NanoStringGeoMxSet object to summarize |
functionList |
optional list of additional functions to calculate negative probe stats, list element names should correspond to expected stat column header |
a NanoStringGeoMxSet object with negative probe summary stats appended to sample data
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) demoData <- summarizeNegatives(demoData, functionList=c(mean=mean, min=min, max=max))
datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") demoData <- readRDS(file.path(datadir, "/demoData.rds")) demoData <- summarizeNegatives(demoData, functionList=c(mean=mean, min=min, max=max))
Update GeoMxSet object to current version
updateGeoMxSet(object)
updateGeoMxSet(object)
object |
GeoMxSet object to update |
updated GeoMxSet object
Take an instance of class NanoStringGeoMxSet
and write
NanoString GeoMx Digital Count Conversion (DCC) data.
writeNanoStringGeoMxSet(x, dir = getwd())
writeNanoStringGeoMxSet(x, dir = getwd())
x |
A NanoStringGeoMxSet object. |
dir |
A directory path to save all the DCC files. |
Zhi Yang & Nicole Ortogero
# Data file paths datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") dccFiles <- dir(datadir, pattern=".dcc$", full.names=TRUE) pkc <- unzip(zipfile = file.path(datadir, "/pkcs.zip")) sampleAnnotationFile <- file.path(datadir, "annotations.xlsx") dccFileColumn <- "Sample_ID" dccSet <- readNanoStringGeoMxSet(dccFiles=dccFiles, pkcFiles=pkc, phenoDataFile=sampleAnnotationFile, phenoDataSheet="CW005", phenoDataDccColName=dccFileColumn, protocolDataColNames=c("aoi", "cell_line", "roi_rep", "pool_rep", "slide_rep"), experimentDataColNames="panel", phenoDataColPrefix="") # All data writeNanoStringGeoMxSet(dccSet)
# Data file paths datadir <- system.file("extdata", "DSP_NGS_Example_Data", package="GeomxTools") dccFiles <- dir(datadir, pattern=".dcc$", full.names=TRUE) pkc <- unzip(zipfile = file.path(datadir, "/pkcs.zip")) sampleAnnotationFile <- file.path(datadir, "annotations.xlsx") dccFileColumn <- "Sample_ID" dccSet <- readNanoStringGeoMxSet(dccFiles=dccFiles, pkcFiles=pkc, phenoDataFile=sampleAnnotationFile, phenoDataSheet="CW005", phenoDataDccColName=dccFileColumn, protocolDataColNames=c("aoi", "cell_line", "roi_rep", "pool_rep", "slide_rep"), experimentDataColNames="panel", phenoDataColPrefix="") # All data writeNanoStringGeoMxSet(dccSet)