TENxFileList: TENxFileList: Represent groups of files from 10X Genomic

View source: R/TENxFileList-class.R

TENxFileListR Documentation

TENxFileList: Represent groups of files from 10X Genomic

Description

This constructor function is meant to handle .tar.gz tarball files from 10X Genomics.

Usage

TENxFileList(..., version, compressed = FALSE)

Arguments

...

Typically, a file path to a tarball archive. Can be named arguments corresponding to file paths, or a named list of file paths.

version

character(1) The version in the tarball. See details.

compressed

logical(1) Whether or not the file provided is compressed, usually as tar.gz (default FALSE)

Details

These tarballs usually contain three files:

  1. matrix.mtx.gz - the counts matrix

  2. features.tsv.gz - row metadata usually represented as rowData

  3. barcodes.tsv.gz - column names corresponding to cell barcode identifiers If all above files are in the tarball, the import method will provide a SingleCellExperiment. Otherwise, a simple list of imported data is given. Note that version "3" uses 'features.tsv.gz' and version "2" uses 'genes.tsv.gz'. If known, indicate the version argument in the TENxFileList constructor function.

Value

Either a SingleCellExperiment or a list of imported data

Examples


fl <- system.file(
    "extdata", "pbmc_granulocyte_sorted_3k_ff_bc_ex_matrix.tar.gz",
    package = "TENxIO", mustWork = TRUE
)

## Method 1 (tarball)
TENxFileList(fl)

## import() method
import(TENxFileList(fl))

## untar to simulate folder output
dir.create(tdir <- tempfile())
untar(fl, exdir = tdir)

## Method 2 (folder)
TENxFileList(tdir)
import(TENxFileList(tdir))

## Method 3 (list of TENxFile objects)
files <- list.files(tdir, recursive = TRUE, full.names = TRUE)
names(files) <- basename(files)
filelist <- lapply(files, TENxFile)

TENxFileList(filelist, compressed = FALSE)

## Method 4 (SimpleList)
TENxFileList(as(filelist, "SimpleList"), compressed = FALSE)

## Method 5 (named arguments)
TENxFileList(
    barcodes.tsv.gz = TENxFile(files[1]),
    features.tsv.gz = TENxFile(files[2]),
    matrix.mtx.gz = TENxFile(files[3])
)

unlink(tdir, recursive = TRUE)


LiNk-NY/TENxIO documentation built on May 3, 2024, 11:08 p.m.