bcbioSingleCell: bcbio single-cell RNA-seq data set

View source: R/AllGenerators.R

bcbioSingleCellR Documentation

bcbio single-cell RNA-seq data set

Description

bcbioSingleCell is an S4 class that extends SingleCellExperiment, and is designed to store a bcbio single-cell RNA-seq analysis. This class contains read counts saved as a sparse matrix (sparseMatrix), sample metadata, and cell quality control metrics.

Usage

bcbioSingleCell(
  uploadDir,
  sampleMetadataFile = NULL,
  organism = NULL,
  ensemblRelease = NULL,
  genomeBuild = NULL,
  gffFile = NULL,
  transgeneNames = NULL,
  interestingGroups = "sampleName"
)

Arguments

uploadDir

character(1). Final upload directory path.

sampleMetadataFile

character(1). Sample metadata file path. CSV or TSV is preferred, but Excel worksheets are also supported. Check the documentation for conventions and required columns.

organism

character(1). Full Latin organism name (e.g. "Homo sapiens").

ensemblRelease

integer(1). Ensembl release version (e.g. 100). We recommend setting this value if possible, for improved reproducibility. When left unset, the latest release available via AnnotationHub/ensembldb is used. Note that the latest version available can vary, depending on the versions of AnnotationHub and ensembldb in use.

genomeBuild

character(1). Ensembl genome build assembly name (e.g. "GRCh38"). If set NULL, defaults to the most recent build available. Note: don't pass in UCSC build IDs (e.g. "hg38").

gffFile

character(1). GFF/GTF (General Feature Format) file. Generally, we recommend using GTF (GFFv2) instead of GFFv3.

transgeneNames

character. Vector indicating which assay rows denote transgenes (e.g. EGFP, TDTOMATO).

interestingGroups

character. Groups of interest to use for visualization. Corresponds to factors describing the columns of the object.

Value

bcbioSingleCell.

Remote data

When working in RStudio, we recommend connecting to the bcbio-nextgen run directory as a remote connection over sshfs.

Note

Updated 2023-09-21.

Author(s)

Michael Steinbaugh

See Also

  • SingleCellExperiment::SingleCellExperiment().

  • .S4methods(class = "bcbioSingleCell").

Examples

uploadDir <- system.file("extdata/indrops", package = "bcbioSingleCell")

x <- bcbioSingleCell(uploadDir)
print(x)

x <- bcbioSingleCell(
    uploadDir = uploadDir,
    sampleMetadataFile = file.path(uploadDir, "metadata.csv")
)
print(x)

hbc/bcbioSingleCell documentation built on Jan. 13, 2024, 2:35 a.m.