Skip to content

Commit 538835a

Browse files
committed
feat: add bambu.singlecell and bambu.spatial wrapper functions
1 parent cee51d6 commit 538835a

2 files changed

Lines changed: 77 additions & 1 deletion

File tree

NAMESPACE

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,8 @@
11
# Generated by roxygen2: do not edit by hand
22

33
export(bambu)
4+
export(bambu.singlecell)
5+
export(bambu.spatial)
46
export(plotBambu)
57
export(prepareAnnotations)
68
export(readFromGTF)

R/bambu.R

Lines changed: 75 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -317,4 +317,78 @@ bambu <- function(reads, annotations = NULL, genome = NULL, NDR = NULL,
317317
}
318318
return(countsSe)
319319
}
320-
}
320+
}
321+
322+
#' Single-cell wrapper for bambu
323+
#' @title Single-cell isoform reconstruction and quantification
324+
#' @description Wrapper around \code{\link{bambu}} for single-cell data.
325+
#' Accepts an \code{opt.singlecell} list for single-cell specific parameters
326+
#' and passes all other arguments directly to \code{\link{bambu}}.
327+
#' @param opt.singlecell A list of single-cell specific parameters:
328+
#' \describe{
329+
#' \item{extractBarcodeUMI}{Logical, whether to extract cell barcodes and
330+
#' UMIs from BAM tags or read names. Defaults to FALSE}
331+
#' \item{dedupUMI}{Logical, whether to perform UMI-based deduplication per
332+
#' barcode. Defaults to FALSE}
333+
#' \item{clusters}{A named list mapping cluster names to barcode vectors,
334+
#' used for cluster-level transcript discovery. Defaults to NULL}
335+
#' }
336+
#' @inheritParams bambu
337+
#' @export
338+
bambu.singlecell <- function(reads, annotations = NULL, genome = NULL, NDR = NULL,
339+
mode = NULL, opt.discovery = NULL, opt.em = NULL, opt.singlecell = NULL,
340+
rcOutDir = NULL, discovery = TRUE, assignDist = TRUE, quant = TRUE,
341+
stranded = FALSE, ncore = 1, yieldSize = NULL, trackReads = FALSE,
342+
returnDistTable = FALSE, lowMemory = FALSE, sampleData = NULL,
343+
fusionMode = FALSE, verbose = FALSE, quantData = NULL,
344+
processByChromosome = FALSE, processByBam = TRUE) {
345+
extractBarcodeUMI <- isTRUE(opt.singlecell$extractBarcodeUMI)
346+
dedupUMI <- isTRUE(opt.singlecell$dedupUMI)
347+
clusters <- opt.singlecell$clusters
348+
bambu(reads = reads, annotations = annotations, genome = genome, NDR = NDR,
349+
mode = mode, opt.discovery = opt.discovery, opt.em = opt.em,
350+
rcOutDir = rcOutDir, discovery = discovery, assignDist = assignDist,
351+
quant = quant, stranded = stranded, ncore = ncore, yieldSize = yieldSize,
352+
trackReads = trackReads, returnDistTable = returnDistTable,
353+
lowMemory = lowMemory, sampleData = sampleData, fusionMode = fusionMode,
354+
verbose = verbose, extractBarcodeUMI = extractBarcodeUMI,
355+
quantData = quantData, dedupUMI = dedupUMI, clusters = clusters,
356+
processByChromosome = processByChromosome, processByBam = processByBam)
357+
}
358+
359+
#' Spatial wrapper for bambu
360+
#' @title Spatial isoform reconstruction and quantification
361+
#' @description Wrapper around \code{\link{bambu}} for spatial transcriptomics
362+
#' data. Accepts an \code{opt.spatial} list for spatial specific parameters
363+
#' and passes all other arguments directly to \code{\link{bambu}}.
364+
#' @param opt.spatial A list of spatial specific parameters:
365+
#' \describe{
366+
#' \item{extractBarcodeUMI}{Logical, whether to extract spot barcodes and
367+
#' UMIs from BAM tags or read names. Defaults to FALSE}
368+
#' \item{dedupUMI}{Logical, whether to perform UMI-based deduplication per
369+
#' barcode. Defaults to FALSE}
370+
#' \item{clusters}{A named list mapping cluster names to barcode vectors,
371+
#' used for cluster-level transcript discovery. Defaults to NULL}
372+
#' }
373+
#' @inheritParams bambu
374+
#' @export
375+
bambu.spatial <- function(reads, annotations = NULL, genome = NULL, NDR = NULL,
376+
mode = NULL, opt.discovery = NULL, opt.em = NULL, opt.spatial = NULL,
377+
rcOutDir = NULL, discovery = TRUE, assignDist = TRUE, quant = TRUE,
378+
stranded = FALSE, ncore = 1, yieldSize = NULL, trackReads = FALSE,
379+
returnDistTable = FALSE, lowMemory = FALSE, sampleData = NULL,
380+
fusionMode = FALSE, verbose = FALSE, quantData = NULL,
381+
processByChromosome = FALSE, processByBam = TRUE) {
382+
extractBarcodeUMI <- isTRUE(opt.spatial$extractBarcodeUMI)
383+
dedupUMI <- isTRUE(opt.spatial$dedupUMI)
384+
clusters <- opt.spatial$clusters
385+
bambu(reads = reads, annotations = annotations, genome = genome, NDR = NDR,
386+
mode = mode, opt.discovery = opt.discovery, opt.em = opt.em,
387+
rcOutDir = rcOutDir, discovery = discovery, assignDist = assignDist,
388+
quant = quant, stranded = stranded, ncore = ncore, yieldSize = yieldSize,
389+
trackReads = trackReads, returnDistTable = returnDistTable,
390+
lowMemory = lowMemory, sampleData = sampleData, fusionMode = fusionMode,
391+
verbose = verbose, extractBarcodeUMI = extractBarcodeUMI,
392+
quantData = quantData, dedupUMI = dedupUMI, clusters = clusters,
393+
processByChromosome = processByChromosome, processByBam = processByBam)
394+
}

0 commit comments

Comments
 (0)