Skip to contents

10x fragment files come in a bed-like format, with columns chr, start, end, cell_id, and pcr_duplicates. Unlike a standard bed format, the format from cellranger has an inclusive end-coordinate, meaning the end coordinate itself is what should be counted as the tagmentation site, rather than offset by 1.

Usage

open_fragments_10x(path, comment = "#", end_inclusive = TRUE)

write_fragments_10x(
  fragments,
  path,
  end_inclusive = TRUE,
  append_5th_column = FALSE
)

Arguments

path

File path (e.g. fragments.tsv or fragments.tsv.gz)

comment

Skip lines at beginning of file which start with comment string

end_inclusive

Whether the end coordinate of the bed is inclusive – i.e. there was an insertion at the end coordinate rather than the base before the end coordinate. This is the 10x default, though it's not quite standard for the bed file format.

fragments

Input fragments object

append_5th_column

Whether to include 5th column of all 0 for compatibility with 10x fragment file outputs (defaults to 4 columns chr,start,end,cell)

Value

10x fragments file object

Details

open_fragments_10x

No disk operations will take place until the fragments are used in a function

write_fragments_10x

Fragments will be written to disk immediately, then returned in a readable object.

Examples

## Download example fragments from pbmc 500 dataset and save in temp directory
data_dir <- file.path(tempdir(), "frags_10x")
dir.create(data_dir, recursive = TRUE, showWarnings = FALSE)
url_base <- "https://cf.10xgenomics.com/samples/cell-atac/2.0.0/atac_pbmc_500_nextgem/"
frags_file <- "atac_pbmc_500_nextgem_fragments.tsv.gz"
atac_raw_url <- paste0(url_base, frags_file)
if (!file.exists(file.path(data_dir, frags_file))) {
 download.file(atac_raw_url, file.path(data_dir, frags_file), mode="wb")
}

#######################################################################
## open_fragments_10x() example
#######################################################################
frags <- open_fragments_10x(
 file.path(data_dir, frags_file)
)
## A Fragments object imported from 10x will not have cell/chromosome 
## information directly known unless written as a BPCells fragment object
frags
#> IterableFragments object of class "ShiftFragments"
#> 
#> Cells: count unknown
#> Chromosomes: count unknown
#> 
#> Queued Operations:
#> 1. Load 10x fragments file from /tmp/RtmpsGFdDm/frags_10x/atac_pbmc_500_nextgem_fragments.tsv.gz
#> 2. Shift start +0bp, end +1bp

frags %>% write_fragments_dir(
 file.path(data_dir, "demo_frags_from_h5"), 
 overwrite = TRUE
)
#> IterableFragments object of class "FragmentsDir"
#> 
#> Cells: 219780 cells with names CCACGTTCAGTAACTC-1, GCGAGAAGTCCACCAG-1 ... AAACGAAGTTCAGAAA-1
#> Chromosomes: 39 chromosomes with names chr1, chr10 ... KI270713.1
#> 
#> Queued Operations:
#> 1. Read compressed fragments from directory /tmp/RtmpsGFdDm/frags_10x/demo_frags_from_h5


#######################################################################
## write_fragments_10x() example
#######################################################################
frags <- write_fragments_10x(
 frags,
 file.path(data_dir, paste0("new_", frags_file))
)
frags
#> IterableFragments object of class "ShiftFragments"
#> 
#> Cells: count unknown
#> Chromosomes: count unknown
#> 
#> Queued Operations:
#> 1. Load 10x fragments file from /tmp/RtmpsGFdDm/frags_10x/new_atac_pbmc_500_nextgem_fragments.tsv.gz
#> 2. Shift start +0bp, end +1bp