| SplikitObject | R Documentation |
R6 class for splicing analysis in single-cell RNA-seq data. Provides a modern object-oriented interface to splikit functionality while maintaining backward compatibility with existing functions.
The SplikitObject encapsulates the core data structures for splicing analysis:
m1: Inclusion matrix (sparse dgCMatrix)
m2: Exclusion matrix (sparse dgCMatrix)
eventData: Event metadata (data.table)
geneExpression: Optional gene expression matrix
m1Inclusion matrix (dgCMatrix). Rows are events, columns are cells.
m2Exclusion matrix (dgCMatrix). Same dimensions as m1.
eventDataEvent metadata (data.table). One row per event.
geneExpressionOptional gene expression matrix (dgCMatrix).
metadataList containing summary statistics and analysis results.
new()Create a new SplikitObject.
SplikitObject$new( junction_ab = NULL, m1 = NULL, m2 = NULL, eventData = NULL, min_counts = 1, verbose = FALSE )
junction_abA junction abundance object from make_junction_ab().
If provided, m1 and eventData are computed automatically.
m1An existing inclusion matrix (dgCMatrix).
m2An existing exclusion matrix (dgCMatrix).
eventDataA data.table with event metadata.
min_countsMinimum count threshold for filtering events (default: 1).
verbosePrint progress messages (default: FALSE).
A new SplikitObject instance.
makeM2()Compute the M2 exclusion matrix from M1 and eventData.
SplikitObject$makeM2( batch_size = 5000, memory_threshold = 2e+09, force_fast = FALSE, n_threads = 1, use_cpp = TRUE, verbose = FALSE )
batch_sizeNumber of groups per batch for memory management (default: 5000).
memory_thresholdMaximum rows before switching to batched processing.
force_fastForce fast processing regardless of size (default: FALSE).
n_threadsNumber of threads for parallel processing (default: 1).
use_cppUse fast C++ implementation (default: TRUE).
verbosePrint progress messages (default: FALSE).
Self (invisibly), for method chaining.
findVariableEvents()Find highly variable splicing events.
SplikitObject$findVariableEvents( min_row_sum = 50, n_threads = 1, verbose = FALSE )
min_row_sumMinimum row sum threshold for filtering (default: 50).
n_threadsNumber of threads for parallel computation (default: 1).
verbosePrint progress messages (default: FALSE).
A data.table with event names and sum_deviance scores.
findVariableGenes()Find highly variable genes from gene expression data.
SplikitObject$findVariableGenes(method = "vst", n_threads = 1, verbose = FALSE)
methodMethod for variable gene selection: "vst" or "sum_deviance" (default: "vst").
n_threadsNumber of threads for parallel computation (default: 1).
verbosePrint progress messages (default: FALSE).
A data.table with gene names and variability scores.
getPseudoCorrelation()Compute pseudo-correlation between splicing and external data.
SplikitObject$getPseudoCorrelation( ZDB_matrix, metric = "CoxSnell", suppress_warnings = TRUE )
ZDB_matrixDense matrix of external data (e.g., gene expression PCs). Must have same dimensions as m1.
metricR-squared metric: "CoxSnell" or "Nagelkerke" (default: "CoxSnell").
suppress_warningsSuppress computation warnings (default: TRUE).
A data.table with event names, pseudo_correlation, and null_distribution.
subset()Subset the object by events and/or cells.
SplikitObject$subset(events = NULL, cells = NULL)
eventsEvent indices or names to keep.
cellsCell indices or names to keep.
Self (invisibly), for method chaining.
setGeneExpression()Set the gene expression matrix.
SplikitObject$setGeneExpression(gene_matrix)
gene_matrixA gene expression matrix (will be converted to dgCMatrix).
Self (invisibly), for method chaining.
annotateEvents()Annotate events with gene information from a GTF file.
SplikitObject$annotateEvents(GTF_file)
GTF_filePath to a GTF annotation file.
Self (invisibly), for method chaining.
summary()Get a summary of the object.
SplikitObject$summary()
A list with object statistics.
print()Print a human-readable summary of the object.
SplikitObject$print()
deepCopy()Create a deep copy of the object.
SplikitObject$deepCopy()
deepIf TRUE, creates a deep copy of all data.
A new SplikitObject with copied data. Validate input matrices and eventData Ensure M2 is computed Ensure matrix is sparse dgCMatrix Standardized error reporting
clone()The objects of this class are cloneable with this method.
SplikitObject$clone(deep = FALSE)
deepWhether to make a deep clone.
# Create from existing matrices using the toy dataset
toy <- load_toy_M1_M2_object()
obj <- SplikitObject$new(m1 = toy$m1, m2 = toy$m2,
eventData = toy$eventdata)
# Find variable events
hve <- obj$findVariableEvents(min_row_sum = 50)
# Or build M2 from M1 + eventData and chain into feature selection
obj2 <- SplikitObject$new(m1 = toy$m1, eventData = toy$eventdata)
results <- obj2$makeM2()$findVariableEvents()
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.