The single constructor for MatisseObject. Combines a
Seurat object with isoform-resolved splicing data and also
computes PSI in one step (unless defer_psi = TRUE). The operating
mode is detected automatically from the inputs you supply:
Usage
CreateMatisseObject(
seurat,
junction_counts = NULL,
transcript_counts = NULL,
events = NULL,
min_coverage = 5L,
defer_psi = FALSE,
verbose = TRUE
)Arguments
- seurat
A
Seuratobject. Required.- junction_counts
A sparse matrix (dgCMatrix, cells x junctions) of raw per-junction read counts. Row names must match
colnames(seurat). Junction IDs (column names) ideally encode genomic coordinates (chr-start-end-strand,chr:start-end:strand, orchr_start_end_strand); Matisse parses these for sashimi plots. Triggers junction mode. Default:NULL.- transcript_counts
A matrix or sparse matrix (transcripts x cells) of raw transcript-level counts. Stored as
Assay5("isoform")in the Seurat object. Column names must overlap withcolnames(seurat). Triggers transcript mode. Default:NULL.- events
Splice-event annotation. Either:
a character vector of paths to SUPPA2
.ioefile(s) (parsed internally), ora pre-built
data.framewith columnsevent_id,gene_id,chr,strand,event_type,inclusion_features,exclusion_features. The*_featurescolumns hold the IDs that support each role — junction IDs in junction mode, transcript IDs in transcript mode.
Required.
- min_coverage
Integer. Minimum total reads per cell per event for the PSI calculation step. Default:
5L. Forwarded toCalculatePSIwhendefer_psi = FALSE.- defer_psi
Logical. Skip the PSI calculation step at construction. The returned object will have raw counts only (no
"psi"assay). CallCalculatePSIlater. Default:FALSE.- verbose
Logical. Print construction progress. Default:
TRUE.
Details
Junction mode (short-read): pass
junction_countsplus splice events viaevents. Junction counts are stored asAssay5("isoform")inside the Seurat object; PSI values land inAssay5("psi").Transcript mode (long-read): pass
transcript_countsplus splice events viaevents(typically a SUPPA2 IOE file). Transcript counts are stored asAssay5("isoform"); PSI values land inAssay5("psi").
Matisse intentionally does not ship a heuristic event caller — bring your own splice-event annotation (SUPPA2 IOE, rMATS table, or equivalent).
Examples
if (FALSE) { # \dontrun{
library(Seurat)
counts <- matrix(rpois(200, 5), nrow = 20,
dimnames = list(paste0("Gene", 1:20),
paste0("Cell", 1:10)))
seu <- CreateSeuratObject(counts)
# Junction mode (short-read): events as a data.frame
jxn <- make_junction_counts()
obj <- CreateMatisseObject(seurat = seu, junction_counts = jxn,
events = my_event_df)
# Transcript mode (long-read): events as SUPPA2 IOE file paths
tx <- make_transcript_counts()
obj <- CreateMatisseObject(seurat = seu, transcript_counts = tx,
events = "path/to/events.ioe")
} # }