Given an input of genomic regions, and a project, this function builds a RangedSummarizedExperiment object. If the user supplies sample_ids, then the RSE will be subsetted to those ids.

libd_rse(regions, project, sample_ids = NULL, chunksize = 1000,
  bpparam = NULL, verbose = TRUE, verboseLoad = verbose,
  round = FALSE)

Arguments

regions

Input of either a GRanges object or a BED file containing Genomic Regions the user wants to build a RSE object for.

project

A character vector of length 1. Valid options are: bsp1 or bsp2.

sample_ids

By default the full phenotype table is returned. The user can supply a character vector of sample ids (typically LIBD RNA numbers such as R10424, R12195).

chunksize

A single integer vector defining the chunksize to use for computing the coverage matrix. Regions will be split into different chunks which can be useful when using a parallel instance as defined by bpparam.

bpparam

A BiocParallelParam-class instance which will be used to calculate the coverage matrix in parallel. By default, SerialParam-class will be used.

verbose

Provides much more information about the flow of information between the client and server.

verboseLoad

If TRUE basic status updates for loading the data will be printed.

round

If TRUE, the counts are rounded to integers. Set to TRUE if you want to match the defaults of recount::scale_counts().

Value

A RangedSummarizedExperiment-class object with the counts scaled to a library size of 40 million reads. If the project data is stranded, then each sample will have two columns (one for each strand: forward and reverse). To collapse the expression into the sense and antisense directions, use sense_antisense on the resulting object from libd_rse().

Details

This function also makes calls to other functions in the libdRSE package including find_samples to obtain the phenotype table for the project and optional sample_ids input, as well as check_region to check if all of the input regions are valid.

Examples

## A test region that overlaps gene SNAP25 snap25 <- GenomicRanges::GRanges('chr20:10286777-10288069:+') ## Compute the RSE for BrainSEQ Phase II snap25_rse <- libd_rse(regions = snap25, 'bsp2', sample_ids = c('R10424', 'R12195', 'R12381', 'R12259'))
#> 2019-08-16 12:38:06 The 'regions' selected are valid regions. #> Returning regions as a GRanges object.
#> 2019-08-16 12:38:06 processing chromosome chr20
#> 2019-08-16 12:38:07 railMatrix: processing regions 1 to 1
#> 2019-08-16 12:38:07 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R10424_C4KC9ACXX_TAATGCGC_L00.Forward.bw
#> 2019-08-16 12:38:08 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R10424_C4KC9ACXX_TAATGCGC_L00.Reverse.bw
#> 2019-08-16 12:38:08 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12195_C4E9TACXX_CGGCTATG_L00.Forward.bw
#> 2019-08-16 12:38:09 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12195_C4E9TACXX_CGGCTATG_L00.Reverse.bw
#> 2019-08-16 12:38:09 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12381_C4UN4ACXX_TCTCGCGC_L00.Forward.bw
#> 2019-08-16 12:38:10 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12381_C4UN4ACXX_TCTCGCGC_L00.Reverse.bw
#> 2019-08-16 12:38:10 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12259_C4L03ACXX_CGCTCATT_L00.Forward.bw
#> 2019-08-16 12:38:10 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12259_C4L03ACXX_CGCTCATT_L00.Reverse.bw
## Explore the resulting object snap25_rse
#> class: RangedSummarizedExperiment #> dim: 1 8 #> metadata(0): #> assays(1): counts #> rownames: NULL #> rowData names(0): #> colnames(8): R10424.Forward R10424.Reverse ... R12381.Forward #> R12381.Reverse #> colData names(79): SAMPLE_ID FQCbasicStats ... #> analysis_eqtl_interaction bigwig_path
## And then the RSE for BrainSEQ Phase I snap25_rse_p1 <- libd_rse(regions = snap25, 'bsp1', sample_ids = c('R2799', 'R2809'))
#> 2019-08-16 12:38:11 The 'regions' selected are valid regions. #> Returning regions as a GRanges object.
#> 2019-08-16 12:38:12 processing chromosome chr20
#> 2019-08-16 12:38:12 railMatrix: processing regions 1 to 1
#> 2019-08-16 12:38:12 railMatrix: processing file /dcl01/ajaffe/data/lab/brainseq_phase1/preprocessed_data/Coverage/R2799_C0JYLACXX.bw
#> 2019-08-16 12:38:13 railMatrix: processing file /dcl01/ajaffe/data/lab/brainseq_phase1/preprocessed_data/Coverage/R2809_D2A01ACXX.bw
## Scenario For more than one chromosome two_chrs <- GenomicRanges::GRanges(c('chr1:1-1000:+', 'chr2:1-1000:-')) two_chrs_rse <- libd_rse(regions = two_chrs, 'bsp2', sample_ids = c('R10424', 'R12195'))
#> 2019-08-16 12:38:13 The 'regions' selected are valid regions. #> Returning regions as a GRanges object.
#> 2019-08-16 12:38:13 processing chromosome chr1
#> 2019-08-16 12:38:14 railMatrix: processing regions 1 to 1
#> 2019-08-16 12:38:14 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R10424_C4KC9ACXX_TAATGCGC_L00.Forward.bw
#> 2019-08-16 12:38:14 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R10424_C4KC9ACXX_TAATGCGC_L00.Reverse.bw
#> 2019-08-16 12:38:14 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12195_C4E9TACXX_CGGCTATG_L00.Forward.bw
#> 2019-08-16 12:38:15 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12195_C4E9TACXX_CGGCTATG_L00.Reverse.bw
#> 2019-08-16 12:38:15 processing chromosome chr2
#> 2019-08-16 12:38:15 railMatrix: processing regions 1 to 1
#> 2019-08-16 12:38:15 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R10424_C4KC9ACXX_TAATGCGC_L00.Forward.bw
#> 2019-08-16 12:38:15 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R10424_C4KC9ACXX_TAATGCGC_L00.Reverse.bw
#> 2019-08-16 12:38:16 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12195_C4E9TACXX_CGGCTATG_L00.Forward.bw
#> 2019-08-16 12:38:16 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12195_C4E9TACXX_CGGCTATG_L00.Reverse.bw
## Two samples with multiple sequence lanes and one sample that is single ## lane snap25_rse_mult <- libd_rse(regions = snap25, 'bsp2', sample_ids = c('R10424', 'R12258', 'R12283'))
#> 2019-08-16 12:38:16 The 'regions' selected are valid regions. #> Returning regions as a GRanges object.
#> 2019-08-16 12:38:16 processing chromosome chr20
#> 2019-08-16 12:38:16 railMatrix: processing regions 1 to 1
#> 2019-08-16 12:38:16 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R10424_C4KC9ACXX_TAATGCGC_L00.Forward.bw
#> 2019-08-16 12:38:17 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R10424_C4KC9ACXX_TAATGCGC_L00.Reverse.bw
#> 2019-08-16 12:38:17 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12258_C4L03ACXX_TCCGGAGA_L00.Forward.bw
#> 2019-08-16 12:38:17 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12258_C4L03ACXX_TCCGGAGA_L00.Reverse.bw
#> 2019-08-16 12:38:18 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12258_C6185ACXX_TCCGGAGA_L005.Forward.bw
#> 2019-08-16 12:38:18 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12258_C6185ACXX_TCCGGAGA_L005.Reverse.bw
#> 2019-08-16 12:38:19 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12283_C5GKYACXX_AGCGATAG_L00.Forward.bw
#> 2019-08-16 12:38:19 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12283_C5GKYACXX_AGCGATAG_L00.Reverse.bw
#> 2019-08-16 12:38:20 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12283_C603PACXX_AGCGATAG_L001.Forward.bw
#> 2019-08-16 12:38:20 railMatrix: processing file /dcl01/lieber/ajaffe/lab/brainseq_phase2/preprocessed_data/DLPFC_RiboZero/Coverage/R12283_C603PACXX_AGCGATAG_L001.Reverse.bw