SingleCellAssay are a generic container for such data and are simple wrappers around SummarizedExperiment objects. Subclasses exist that embue the container with additional attributes, eg FluidigmAssay.

FromFlatDF(
  dataframe,
  idvars,
  primerid,
  measurement,
  id = numeric(0),
  cellvars = NULL,
  featurevars = NULL,
  phenovars = NULL,
  class = "SingleCellAssay",
  check_sanity = TRUE,
  ...
)

Arguments

dataframe

A 'flattened' data.frame or data.table containing columns giving cell and feature identifiers and a measurement column

idvars

character vector naming columns that uniquely identify a cell

primerid

character vector of length 1 that names the column that identifies what feature (i.e. gene) was measured

measurement

character vector of length 1 that names the column containing the measurement

id

An identifier (eg, experiment name) for the resulting object

cellvars

Character vector naming columns containing additional cellular metadata

featurevars

Character vector naming columns containing additional feature metadata

phenovars

Character vector naming columns containing additional phenotype metadata

class

desired subclass of object. Default SingleCellAssay.

check_sanity

(default: TRUE) Set FALSE to override sanity checks that try to ensure that the default assay is log-transformed and has at least one exact zero. See defaultAssay for details on the "default assay" which is assumed to contain log transformed data.

...

additional arguments are ignored

Value

SingleCellAssay, or derived, object

Examples

data(vbeta)
colnames(vbeta)
#>  [1] "Sample.ID"         "Subject.ID"        "Experiment.Number"
#>  [4] "Chip.Number"       "Stim.Condition"    "Time"             
#>  [7] "Population"        "Number.of.Cells"   "Well"             
#> [10] "Gene"              "Ct"               
vbeta <- computeEtFromCt(vbeta)
vbeta.fa <- FromFlatDF(vbeta, idvars=c("Subject.ID", "Chip.Number", "Well"),
primerid='Gene', measurement='Et', ncells='Number.of.Cells',
geneid="Gene",cellvars=c('Number.of.Cells', 'Population'),
phenovars=c('Stim.Condition','Time'), id='vbeta all', class='FluidigmAssay')
#> No dimnames in `exprsArray`, assuming `fData` and `cData` are sorted according to `exprsArray`
#> Assuming data assay in position 1, with name Et is log-transformed.
show(vbeta.fa)
#> class: FluidigmAssay 
#> dim: 75 456 
#> metadata(0):
#> assays(1): Et
#> rownames(75): B3GAT1 BAX ... TNFRSF9 TNFSF10
#> rowData names(2): Gene primerid
#> colnames(456): Sub01 1 A01 Sub01 1 A02 ... Sub02 3 H10 Sub02 3 H11
#> colData names(9): Number.of.Cells Population ... Time wellKey
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):
nrow(vbeta.fa)
#> [1] 75
ncol(vbeta.fa)
#> [1] 456
head(mcols(vbeta.fa)$primerid)
#> [1] "B3GAT1" "BAX"    "BCL2"   "CCL2"   "CCL3"   "CCL4"  
table(colData(vbeta.fa)$Subject.ID)
#> 
#> Sub01 Sub02 
#>   177   279 
vbeta.sub <- subset(vbeta.fa, Subject.ID=='Sub01')
show(vbeta.sub)
#> class: FluidigmAssay 
#> dim: 75 177 
#> metadata(0):
#> assays(1): Et
#> rownames(75): B3GAT1 BAX ... TNFRSF9 TNFSF10
#> rowData names(2): Gene primerid
#> colnames(177): Sub01 1 A01 Sub01 1 A02 ... Sub01 2 H09 Sub01 2 H10
#> colData names(9): Number.of.Cells Population ... Time wellKey
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(0):