Function to import and convert tabular data from a spreadsheet or
a data.frame
into a SingleCellExperiment
and QFeatures
object.
readSCP(...)
readSCPfromDIANN(...)
readSingleCellExperiment(...)
An instance of class SingleCellExperiment
or a
QFeatures
, composed of SingleCellExperiment
objects.
The SingleCellExperiment
class is built on top of the
RangedSummarizedExperiment
class. This means that some column names
are forbidden in the rowData
. Avoid using the following names:
seqnames
, ranges
, strand
, start
, end
,
width
, element
The more general QFeatures::readQFeatures()
function, which
this function depends on.
The more general QFeatures::readQFeaturesFromDIANN()
function,
for details and an example on how to read label-free and plexDIA
(mTRAQ) data processed with DIA-NN.
The QFeatures::readSummarizedExperiment()
function, which
readSingleCellExperiment()
depends on.
The SingleCellExperiment::SingleCellExperiment()
class.
######################################################
## Load a single acquisition as a SingleCellExperiment
## Load a data.frame with PSM-level data
data("mqScpData")
## Create the QFeatures object
sce <- readSingleCellExperiment(mqScpData,
quantCols = grep("RI", colnames(mqScpData)))
sce
#> class: SingleCellExperiment
#> dim: 1361 0
#> metadata(0):
#> assays(1): ''
#> rownames(1361): 1 2 ... 1360 1361
#> rowData names(0):
#> colnames: NULL
#> colData names(0):
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
######################################################
## Load multiple acquisitions as a QFeatures
## Load an example table containing MaxQuant output
data("mqScpData")
## Load the (user-generated) annotation table
data("sampleAnnotation")
## Format the tables into a QFeatures object
readSCP(assayData = mqScpData,
colData = sampleAnnotation,
runCol = "Raw.file")
#> Checking arguments.
#> Loading data as a 'SummarizedExperiment' object.
#> Splitting data in runs.
#> Formatting sample annotations (colData).
#> Formatting data as a 'QFeatures' object.
#> An instance of class QFeatures containing 4 assays:
#> [1] 190222S_LCA9_X_FP94BM: SingleCellExperiment with 395 rows and 16 columns
#> [2] 190321S_LCA10_X_FP97AG: SingleCellExperiment with 487 rows and 16 columns
#> [3] 190321S_LCA10_X_FP97_blank_01: SingleCellExperiment with 109 rows and 16 columns
#> [4] 190914S_LCB3_X_16plex_Set_21: SingleCellExperiment with 370 rows and 16 columns