Using output from fullCoverage, export the coverage from all the samples to BigWig files using createBwSample.

createBw(fullCov, path = ".", keepGR = TRUE, ...)

Arguments

fullCov

A list where each element is the result from loadCoverage used with returnCoverage = TRUE. Can be generated using fullCoverage.

path

The path where the BigWig files will be created.

keepGR

If TRUE, the GRanges objects created by coerceGR grouped into a GRangesList are returned. Otherwise they are discarded.

...

Arguments passed to createBwSample.

Value

If keepGR = TRUE, then a GRangesList

with the output for coerceGR for each of the samples.

Details

Use at most one core per chromosome.

Author

Leonardo Collado-Torres

Examples

## Create a small fullCov object with data only for chr21
fullCov <- list("chr21" = genomeDataRaw)

## Keep only 2 samples
fullCov$chr21$coverage <- fullCov$chr21$coverage[c(1, 31)]

## Create the BigWig files for all samples in a test dir
dir.create("createBw-example")
bws <- createBw(fullCov, "createBw-example")
#> 2023-05-07 06:01:15.845656 coerceGR: coercing sample ERR009101
#> 2023-05-07 06:01:15.876899 createBwSample: exporting bw for sample ERR009101
#> 2023-05-07 06:01:15.915947 coerceGR: coercing sample SRR031960
#> 2023-05-07 06:01:15.94209 createBwSample: exporting bw for sample SRR031960
#> Warning: There are no bases with coverage > 0 for sample SRR031960. Thus no bw file will be created for this sample.

## Explore the output
bws
#> GRangesList object of length 2:
#> $ERR009101
#> GRanges object with 16 ranges and 1 metadata column:
#>         seqnames            ranges strand |     score
#>            <Rle>         <IRanges>  <Rle> | <integer>
#>   chr21    chr21 47410303-47410336      * |         1
#>   chr21    chr21 47410687-47410688      * |         1
#>   chr21    chr21 47410699-47410734      * |         1
#>   chr21    chr21 47410951-47410955      * |         1
#>   chr21    chr21          47411924      * |         1
#>     ...      ...               ...    ... .       ...
#>   chr21    chr21 47415237-47415272      * |         1
#>   chr21    chr21 47415317-47415352      * |         1
#>   chr21    chr21 47417629-47417664      * |         1
#>   chr21    chr21 47417674-47417676      * |         1
#>   chr21    chr21 47418035-47418067      * |         1
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome
#> 
#> $SRR031960
#> GRanges object with 0 ranges and 1 metadata column:
#>    seqnames    ranges strand |     score
#>       <Rle> <IRanges>  <Rle> | <integer>
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome
#> 

## First sample
bws[[1]]
#> GRanges object with 16 ranges and 1 metadata column:
#>         seqnames            ranges strand |     score
#>            <Rle>         <IRanges>  <Rle> | <integer>
#>   chr21    chr21 47410303-47410336      * |         1
#>   chr21    chr21 47410687-47410688      * |         1
#>   chr21    chr21 47410699-47410734      * |         1
#>   chr21    chr21 47410951-47410955      * |         1
#>   chr21    chr21          47411924      * |         1
#>     ...      ...               ...    ... .       ...
#>   chr21    chr21 47415237-47415272      * |         1
#>   chr21    chr21 47415317-47415352      * |         1
#>   chr21    chr21 47417629-47417664      * |         1
#>   chr21    chr21 47417674-47417676      * |         1
#>   chr21    chr21 47418035-47418067      * |         1
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome

## Note that if a sample has no bases with coverage > 0, the GRanges object
## is empty and no BigWig file is created for that sample.
bws[[2]]
#> GRanges object with 0 ranges and 1 metadata column:
#>    seqnames    ranges strand |     score
#>       <Rle> <IRanges>  <Rle> | <integer>
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome


## Exports fullCoverage() output to BigWig files