extract_design_mdesign_sbs96: Extract Design-Metadesign (XU) matrix for SBS-96 meta-feature...

View source: R/extract_design_mdesign_sbs96.R

extract_design_mdesign_sbs96R Documentation

Extract Design-Metadesign (XU) matrix for SBS-96 meta-feature categories from a maf file

Description

Extract Design-Metadesign (XU) matrix for SBS-96 meta-feature categories from a maf file

Usage

extract_design_mdesign_sbs96(
  maf,
  chromosome_col = "Chromosome",
  start_position_col = "Start_Position",
  end_position_col = "End_Position",
  ref_col = "Reference_Allele",
  alt_col = "Tumor_Seq_Allele2",
  sample_id_col = "sample",
  return_vranges_obj = FALSE,
  ...
)

Arguments

maf

mutation annotation file – a data frame-like object with at least six columns for: Chromosome, Start_position, End_position, Reference_Allele, Tumor_Seq_Allele2 and Sample ID. NOTE: uniqueness of rows of maf is assumed.

sample_id_col

name of the column in maf containing tumor sample IDs.

...

Unused.

return_vranges_object

logical. Should the intermediate 'VRanges' object from {SomaticSignatures} cataloging SBS-96 contexts for all mutations in maf be returned as an attribute named VRanges of the output? Defaults to FALSE.

Details

This function calls uses functions VariantAnnotation::VRanges(), IRanges::IRanges(), SomaticSignatures::ucsc() and SomaticSignatures::mutationContext(), together with the BSgenome.Hsapiens.UCSC.hg19 database to obtain the 96 single base substitution contexts for each mutation in the maf file in the form of a 'VRanges' object from SomaticSignatures. Then using SomaticSignatures::motifMatrix a 96 x n_tumors is calculated, which is subsequently transposed and converted into a sparseMatrix object in the form of a n_patient x 96 design matrix

Value

An n_tumor x 96 sparse dgCMatrix, with (i, j)th entry providing the total number of variants in tumor i associated with j-th SBS-96 category.

Note

The bioconductor packages SomaticSignatures and VariantAnnotation are not automatically installed with hidgenclassifier. Please install them separately.

Using any SomaticSignatures function triggers the loading of proxy and GGally packages, which overwrites defaults S3 methods of a few functions from ggplot2 and registry.

Author(s)

Ronglai Shen, Saptarshi Chakraborty

Examples

data("impact")
sbs96_mdesign <- extract_design_mdesign_sbs96(
  maf = impact,
  chromosome_col = "Chromosome",
  start_position_col = "Start_Position",
  end_position_col = "End_Position",
  ref_col = "Reference_Allele",
  alt_col = "Tumor_Seq_Allele2",
  sample_id_col = "patient_id"
)
dim(sbs96_mdesign)


c7rishi/hidgenclassifier documentation built on June 14, 2024, 11:10 a.m.