Skip to contents

makeGRangesListFromCopyNumber allows the user to convert objects of class data.frame or S4Vectors::DataFrame to a GRangesList. It includes additional features specific to TCGA data such as, hugo symbols, probe numbers, segment means, and ucsc build (if available).

Usage

makeGRangesListFromCopyNumber(
  df,
  split.field,
  names.field = "Hugo_Symbol",
  ...
)

Arguments

df

A data.frame or DataFrame class object. list class objects are coerced to data.frame or DataFrame.

split.field

A character vector of length one indicating the column to be used as sample identifiers

names.field

A character vector of length one indicating the column to be used as names for each of the ranges in the data

...

Additional arguments to pass on to GenomicRanges::makeGRangesListFromDataFrame

Value

A GRangesList class object

Examples

library(GenomicDataCommons)

manif <- files() |>
    filter(~ cases.project.project_id == "TCGA-COAD" &
        data_type == "Copy Number Segment") |>
    manifest(size = 1)

fname <- gdcdata(manif$id)

UUIDtoBarcode(names(fname), from_type = "file_id")
#>                                file_id associated_entities.entity_submitter_id
#> 1 7221548d-4204-4fdb-b82c-3224cc15518d            TCGA-AA-3715-10B-01D-A91V-36
#> 2 7221548d-4204-4fdb-b82c-3224cc15518d            TCGA-AA-3715-01A-01D-0957-02

cndata <- read.delim(fname[[1L]])

makeGRangesListFromCopyNumber(
    df = cndata,
    split.field = "GDC_Aliquot_ID",
    keep.extra.columns = TRUE
)
#> GRangesList object of length 1:
#> $`TCGA-AA-3715-01A-01D-0957-02`
#> GRanges object with 6167 ranges and 2 metadata columns:
#>          seqnames            ranges strand | Num_Probes Segment_Mean
#>             <Rle>         <IRanges>  <Rle> |  <integer>    <numeric>
#>      [1]     chr1      17001-635000      * |         13    -2.076454
#>      [2]     chr1     701001-771000      * |          3   -29.246122
#>      [3]     chr1     777001-811000      * |         29    -1.710666
#>      [4]     chr1    811001-1649000      * |        800    -0.371095
#>      [5]     chr1   1649001-3687000      * |       1874    -0.254668
#>      ...      ...               ...    ... .        ...          ...
#>   [6163]     chrY 56855001-56863000      * |          8    -0.571704
#>   [6164]     chrY 56863001-56865000      * |          2    -3.723692
#>   [6165]     chrY 56865001-56876000      * |         11    -1.073070
#>   [6166]     chrY 56876001-56887000      * |         11    -0.310591
#>   [6167]     chrM           1-16569      * |         17     0.763246
#>   -------
#>   seqinfo: 25 sequences from an unspecified genome; no seqlengths
#>