SplikitObject: SplikitObject

SplikitObjectR Documentation

SplikitObject

Description

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.

Details

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

Public fields

m1

Inclusion matrix (dgCMatrix). Rows are events, columns are cells.

m2

Exclusion matrix (dgCMatrix). Same dimensions as m1.

eventData

Event metadata (data.table). One row per event.

geneExpression

Optional gene expression matrix (dgCMatrix).

metadata

List containing summary statistics and analysis results.

Methods

Public methods


Method new()

Create a new SplikitObject.

Usage
SplikitObject$new(
  junction_ab = NULL,
  m1 = NULL,
  m2 = NULL,
  eventData = NULL,
  min_counts = 1,
  verbose = FALSE
)
Arguments
junction_ab

A junction abundance object from make_junction_ab(). If provided, m1 and eventData are computed automatically.

m1

An existing inclusion matrix (dgCMatrix).

m2

An existing exclusion matrix (dgCMatrix).

eventData

A data.table with event metadata.

min_counts

Minimum count threshold for filtering events (default: 1).

verbose

Print progress messages (default: FALSE).

Returns

A new SplikitObject instance.


Method makeM2()

Compute the M2 exclusion matrix from M1 and eventData.

Usage
SplikitObject$makeM2(
  batch_size = 5000,
  memory_threshold = 2e+09,
  force_fast = FALSE,
  n_threads = 1,
  use_cpp = TRUE,
  verbose = FALSE
)
Arguments
batch_size

Number of groups per batch for memory management (default: 5000).

memory_threshold

Maximum rows before switching to batched processing.

force_fast

Force fast processing regardless of size (default: FALSE).

n_threads

Number of threads for parallel processing (default: 1).

use_cpp

Use fast C++ implementation (default: TRUE).

verbose

Print progress messages (default: FALSE).

Returns

Self (invisibly), for method chaining.


Method findVariableEvents()

Find highly variable splicing events.

Usage
SplikitObject$findVariableEvents(
  min_row_sum = 50,
  n_threads = 1,
  verbose = FALSE
)
Arguments
min_row_sum

Minimum row sum threshold for filtering (default: 50).

n_threads

Number of threads for parallel computation (default: 1).

verbose

Print progress messages (default: FALSE).

Returns

A data.table with event names and sum_deviance scores.


Method findVariableGenes()

Find highly variable genes from gene expression data.

Usage
SplikitObject$findVariableGenes(method = "vst", n_threads = 1, verbose = FALSE)
Arguments
method

Method for variable gene selection: "vst" or "sum_deviance" (default: "vst").

n_threads

Number of threads for parallel computation (default: 1).

verbose

Print progress messages (default: FALSE).

Returns

A data.table with gene names and variability scores.


Method getPseudoCorrelation()

Compute pseudo-correlation between splicing and external data.

Usage
SplikitObject$getPseudoCorrelation(
  ZDB_matrix,
  metric = "CoxSnell",
  suppress_warnings = TRUE
)
Arguments
ZDB_matrix

Dense matrix of external data (e.g., gene expression PCs). Must have same dimensions as m1.

metric

R-squared metric: "CoxSnell" or "Nagelkerke" (default: "CoxSnell").

suppress_warnings

Suppress computation warnings (default: TRUE).

Returns

A data.table with event names, pseudo_correlation, and null_distribution.


Method subset()

Subset the object by events and/or cells.

Usage
SplikitObject$subset(events = NULL, cells = NULL)
Arguments
events

Event indices or names to keep.

cells

Cell indices or names to keep.

Returns

Self (invisibly), for method chaining.


Method setGeneExpression()

Set the gene expression matrix.

Usage
SplikitObject$setGeneExpression(gene_matrix)
Arguments
gene_matrix

A gene expression matrix (will be converted to dgCMatrix).

Returns

Self (invisibly), for method chaining.


Method annotateEvents()

Annotate events with gene information from a GTF file.

Usage
SplikitObject$annotateEvents(GTF_file)
Arguments
GTF_file

Path to a GTF annotation file.

Returns

Self (invisibly), for method chaining.


Method summary()

Get a summary of the object.

Usage
SplikitObject$summary()
Returns

A list with object statistics.


Method print()

Print a human-readable summary of the object.

Usage
SplikitObject$print()

Method deepCopy()

Create a deep copy of the object.

Usage
SplikitObject$deepCopy()
Arguments
deep

If TRUE, creates a deep copy of all data.

Returns

A new SplikitObject with copied data. Validate input matrices and eventData Ensure M2 is computed Ensure matrix is sparse dgCMatrix Standardized error reporting


Method clone()

The objects of this class are cloneable with this method.

Usage
SplikitObject$clone(deep = FALSE)
Arguments
deep

Whether to make a deep clone.

Examples


# 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()



splikit documentation built on May 13, 2026, 9:08 a.m.