View source: R/read10xVisium.R
read10xVisium | R Documentation |
Creates a SpatialExperiment
from the Space Ranger
output directories for 10x Genomics Visium spatial gene expression data.
read10xVisium(
samples = "",
sample_id = paste0("sample", sprintf("%02d", seq_along(samples))),
type = c("HDF5", "sparse"),
data = c("filtered", "raw"),
images = "lowres",
load = TRUE
)
samples |
a character vector specifying one or more directories, each corresponding to a 10x Genomics Visium sample (see Details); if provided, names will be used as sample identifiers |
sample_id |
character string specifying unique sample identifiers,
one for each directory specified via |
type |
character string specifying
the type of format to read count data from
(see |
data |
character string specifying whether to read in filtered (spots mapped to tissue) or raw data (all spots). |
images |
character vector specifying which images to include.
Valid values are |
load |
logical; should the image(s) be loaded into memory
as a |
The constructor assumes data from each sample are located in a single output directory as returned by Space Ranger, thus having the following file organization (where "raw/filtered" refers to either "raw" or "filtered" to match the 'data' argument.) The base directory "outs/" from Space Ranger can either be included manually in the paths provided in 'samples', or can be ignored; if ignored, it will be added automatically. The '.h5' files are used if 'type = "HDF5"'. (Note that 'tissue_positions.csv' was renamed in Space Ranger v2.0.0.)
sample
· | — outs
· · | — raw/filtered_feature_bc_matrix.h5
· · | — raw/filtered_feature_bc_matrix
· · · | — barcodes.tsv.gz
· · · | — features.tsv.gz
· · · | — matrix.mtx.gz
· · | — spatial
· · · | — scalefactors_json.json
· · · | — tissue_lowres_image.png
· · · | — tissue_positions.csv
a SpatialExperiment
object
Helena L. Crowell
dir <- system.file(
file.path("extdata", "10xVisium"),
package = "SpatialExperiment")
sample_ids <- c("section1", "section2")
samples <- file.path(dir, sample_ids, "outs")
list.files(samples[1])
list.files(file.path(samples[1], "spatial"))
file.path(samples[1], "raw_feature_bc_matrix")
(spe <- read10xVisium(samples, sample_ids,
type = "sparse", data = "raw",
images = "lowres", load = FALSE))
# base directory 'outs/' from Space Ranger can also be omitted
samples2 <- file.path(dir, sample_ids)
(spe2 <- read10xVisium(samples2, sample_ids,
type = "sparse", data = "raw",
images = "lowres", load = FALSE))
# tabulate number of spots mapped to tissue
cd <- colData(spe)
table(
in_tissue = cd$in_tissue,
sample_id = cd$sample_id)
# view available images
imgData(spe)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.