Title: | Gene Break Detection |
---|---|
Description: | Recurrent breakpoint gene detection on copy number aberration profiles. |
Authors: | Evert van den Broek, Stef van Lieshout |
Maintainer: | Evert van den Broek <[email protected]> |
License: | GPL-2 |
Version: | 1.37.0 |
Built: | 2024-10-30 07:22:01 UTC |
Source: | https://github.com/bioc/GeneBreak |
Access Object Data. This method lists possible functions to access the data of the object.
## S4 method for signature 'CopyNumberBreakPoints' accessOptions(object)
## S4 method for signature 'CopyNumberBreakPoints' accessOptions(object)
object |
An object of class |
prints text to screen
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) accessOptions( bp )
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) accessOptions( bp )
Maps features to gene locations.
## S4 method for signature 'CopyNumberBreakPoints' addGeneAnnotation(object, geneAnnotation)
## S4 method for signature 'CopyNumberBreakPoints' addGeneAnnotation(object, geneAnnotation)
object |
An object of class CopyNumberBreakPoints |
geneAnnotation |
An object of class |
The end of the first feature after gene start location up to and including the first feature after gene end location will be defined as gene-associated feaures. For hg18, hg19 and hg38 built-in gene annotation files obtained from ensembl can be used. Please take care to use a matching reference genome for your breakpoint data. In stead of using the built-in gene annotion files, feature-to-gene mapping can be based on an user-defined annotion file. The dataframe should contain at least these four columns: "Gene", "Chromosome", "Start" and "End".
Returns an object of class CopyNumberBreakPointGenes with gene annotation added.
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) ## other buil-in gene annotations are: # data( ens.gene.ann.hg19 ) # data( ens.gene.ann.hg38 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) # input copynumber.data.chr20 is hg18 based bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) ## options to inspect the data bp accessOptions( bp )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) ## other buil-in gene annotations are: # data( ens.gene.ann.hg19 ) # data( ens.gene.ann.hg38 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) # input copynumber.data.chr20 is hg18 based bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) ## options to inspect the data bp accessOptions( bp )
Selects breakpoints by filter criteria options.
## S4 method for signature 'CopyNumberBreakPoints' bpFilter(object, filter = "CNA-ass", threshold = NULL)
## S4 method for signature 'CopyNumberBreakPoints' bpFilter(object, filter = "CNA-ass", threshold = NULL)
object |
An object of class CopyNumberBreakPoints |
filter |
Type of filter. This can be either "CNA-ass", "deltaSeg" or "deltaCall".
|
threshold |
Set the minimal log2 ratio difference between segments. This parameter is required for the "deltaSeg" filter option |
Filter options "CNA-ass" and "deltaCall" require calls in addition to segmented copynumber data (see input for getBreakpoints()
)
Returns an object of class CopyNumberBreakPoints with breakpoint matrix replaced by filtered breakpoints.
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp, filter = "CNA-ass" ) bp <- bpFilter( bp, filter = "deltaSeg", threshold = 0.2 ) ## options to inspect the data bp accessOptions( bp )
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp, filter = "CNA-ass" ) bp <- bpFilter( bp, filter = "deltaSeg", threshold = 0.2 ) ## options to inspect the data bp accessOptions( bp )
Indentifies genes affected by breakpoint locations.
## S4 method for signature 'CopyNumberBreakPointGenes' bpGenes(object)
## S4 method for signature 'CopyNumberBreakPointGenes' bpGenes(object)
object |
An object of class CopyNumberBreakPointGenes |
This step requires feature-to-gene annotations added to the input object (see ?addGeneAnnotation
).
Returns an object of class CopyNumberBreakPointGenes with gene-breakpoint information
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) ## options to inspect the data bp accessOptions( bp )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) ## options to inspect the data bp accessOptions( bp )
Plots breakpoint frequencies per chromosome
## S4 method for signature 'CopyNumberBreakPoints' bpPlot(object, plot.chr = NULL, plot.ylim = 15, fdr.threshold = 0.1, add.jitter = FALSE)
## S4 method for signature 'CopyNumberBreakPoints' bpPlot(object, plot.chr = NULL, plot.ylim = 15, fdr.threshold = 0.1, add.jitter = FALSE)
object |
An object of class CopyNumberBreakPoints or CopyNumberBreakPointGenes |
plot.chr |
A vector with chromosome(s) to plot. All chromosomes will be plotted when NULL is used. |
plot.ylim |
An integer giving the max y coordinate. |
fdr.threshold |
The FDR threshold to label recurrent breakpoint genes with their gene name |
add.jitter |
Logical. If TRUE, function jitter will be used for the y position of gene labels |
The plot includes breakpoint locations and breakpoint gene frequencies. Genes that are recurrently affected are labeled with their gene name.
calls plot function
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) bp <- bpStats( bp ) bpPlot( bp, c(20) )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) bp <- bpStats( bp ) bpPlot( bp, c(20) )
Applies cohort-based statistics to identify genes and/or chromosomal locations that are recurrently affected by breakpoints.
## S4 method for signature 'CopyNumberBreakPoints' bpStats(object, level = "gene", method = "BH", fdr.threshold = 1)
## S4 method for signature 'CopyNumberBreakPoints' bpStats(object, level = "gene", method = "BH", fdr.threshold = 1)
object |
An object of class CopyNumberBreakPointGenes |
level |
The level at which to operate, this can be either "gene" (correcting for gene length) or "feature" (per probe/bin) |
method |
The FDR correction method to apply. This can be "BH" (applies Benjamini-Hochberg-type FDR correction) or "Gilbert" (for dedicated Benjamini-Hochberg-type FDR correction) |
fdr.threshold |
The threshold for FDR correction' |
The statistical method on gene-level corrects for covariates that may influence the probability to be a breakpoint gene including number of breakpoints in a profile, number of gene-associated features and gene length by gene-associated feature coverage. The statistical analysis includes multiple testing where standard Benjamini-Hochberg-type FDR correction will be performed by default. This less computational intensive method assumes a similar null-distribution for all candidate breakpoint events and satisfies for analysis on breakpoint location-level. For statistics on gene-level however, we recommend to apply the more comprehensive and powerful dedicated Benjamini-Hochberg-type FDR correction that accounts for discreteness in null-distribution (Gilbert, 2005) following correction for covariates that may influence the probability to be a breakpoint gene including number of breakpoints in a profile, number of gene-associated features and gene length by gene-associated feature coverage.
Returns an object of class CopyNumberBreakPointGenes with cohort based statistics added.
Gilbert,P.B. (2005) A modified false discovery rate multiple-comparisons procedure for discrete data, applied to human immunodeficiency virus genetics. Journal of the Royal Statistical Society Series C-Applied Statistics, 54, 143-158.
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) bp <- bpStats( bp ) ## options to inspect the data bp accessOptions( bp )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) bp <- bpStats( bp ) ## options to inspect the data bp accessOptions( bp )
Access Object breakpointData. This method returns a dataframe with breakpoint values per feature.
## S4 method for signature 'CopyNumberBreakPoints' breakpointData(object)
## S4 method for signature 'CopyNumberBreakPoints' breakpointData(object)
object |
An object of class |
a dataframe with breakpoint values
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) breakpointData( bp )
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) breakpointData( bp )
Access Object breakpointsPerGene. This method returns a dataframe with breakpoints per gene.
## S4 method for signature 'CopyNumberBreakPointGenes' breakpointsPerGene(object)
## S4 method for signature 'CopyNumberBreakPointGenes' breakpointsPerGene(object)
object |
An object of class |
a dataframe with breakpoints per gene
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) breakpointsPerGene( bp )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) breakpointsPerGene( bp )
Access Object callData. This method returns a dataframe with feature call values.
## S4 method for signature 'CopyNumberBreakPoints' callData(object)
## S4 method for signature 'CopyNumberBreakPoints' callData(object)
object |
An object of class |
a dataframe with feature call values
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) callData( bp )
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) callData( bp )
A test dataset containing copynumber data of chromosome 18 for the GeneBreak package (hg18 based). This copy number aberration (CNA) data was obtained by analysis of 200 array-CGH (Agilent 180k) samples from advanced colorectal cancers.
data( copynumber.data.chr18 )
data( copynumber.data.chr18 )
An object of class cghCall
An object of class cghCall
A test dataset containing chromosome 20 copynumber data for the GeneBreak package (hg18 based). This copy number aberration (CNA) data was obtained by analysis of 200 array-CGH (Agilent 180k) samples from advanced colorectal cancers.
data( copynumber.data.chr20 )
data( copynumber.data.chr20 )
An object of class cghCall
An object of class cghCall
A test dataset containing chromosome 21 copynumber data for the GeneBreak package (hg18 based). This copy number aberration (CNA) data was obtained by analysis of 200 array-CGH (Agilent 180k) samples from advanced colorectal cancers.
data( copynumber.data.chr21 )
data( copynumber.data.chr21 )
An object of class cghCall
An object of class cghCall
An S4 class to represent a CopyNumberBreakPointGenes object
geneAnnotation
A data.frame with original gene annotation input
geneData
A data.frame with gene information added by package methods
featuresPerGene
A list with the associated features per gene
breakpointsPerGene
A matrix with breakage status per gene
callData( object )
Returns feature call values:
segmentData( object )
Returns feature segment values
breakpointData( object )
Returns feature breakpoint values
sampleNames( object )
Returns vector with sample names
namesFeatures( object )
Returns vector with feature names
featureChromosomes( object )
Returns vector of feature chromosomes
featureInfo( object )
Returns feature data/information
geneChromosomes( object )
Returns vector of gene chromosomes
geneInfo( object )
Returns gene data/information
featuresPerGene( object )
Returns a list of genes with coupled features
breakpointsPerGene( object )
Returns gene break status
recurrentGenes( object )
Returns recurrently broken genes
getBreakpoints Builds the CopyNumberBreakPoints object from copynumber data and detects breakpoint locations
bpFilter Selects breakpoints by filter criteria options
addGeneAnnotation Maps features to gene locations
bpGenes Indentifies genes affected by breakpoint locations
bpStats Applies cohort-based statistics to identify genes and/or chromosomal locations that are recurrently affected by breakpoints
bpPlot Plots breakpoint frequencies per chromosome
E. van den Broek and S. van Lieshout
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) bp <- bpStats( bp ) bpPlot( bp, c(20) )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) bp <- bpStats( bp ) bpPlot( bp, c(20) )
An S4 class to represent a CopyNumberBreakPoints object.
segmDiff
A matrix with breakpoints based on segment values
callDiff
A matrix with breakpoints based on call values
segments
A matrix with segmented copy number values
calls
A matrix with copy number calls
featureAnnotation
A dataframe with predefined information about the features (usually probes or bins)
featureData
A dataframe with calculated information about the features (usually probes or bins)
callData( object ) Returns feature call values
segmentData( object ) Returns feature segment values
breakpointData( object ) Returns feature breakpoint values
sampleNames( object ) Returns vector with sample names
namesFeatures( object ) Returns vector with feature names
featureChromosomes( object ) Returns vector of feature chromosomes
featureInfo( object ) Returns feature data/information
getBreakpoints Builds the CopyNumberBreakPoints object from copynumber data and detects breakpoint locations
bpFilter Selects breakpoints by filter criteria options
bpStats Applies cohort-based statistics to identify chromosomal locations that are recurrently affected by breakpoints
bpPlot Plots breakpoint frequencies per chromosome
E. van den Broek and S. van Lieshout
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- bpStats( bp , level = 'feature' , method = 'BH' ) bpPlot( bp, c(20) )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- bpStats( bp , level = 'feature' , method = 'BH' ) bpPlot( bp, c(20) )
A dataset containing the gene locations based on human genome reference hg18 that was obtained from BioMart.
data( ens.gene.ann.hg18 )
data( ens.gene.ann.hg18 )
A data.frame
Dataframe with 5 columns:
Gene: ensembl gene name
EnsID: ensembl gene id
Chromosome: Genomic Chromosome
Start: Genomic start of gene
End: Genomic end of gene
data.frame
A dataset containing the gene locations based on human genome reference hg19 that was obtained from BioMart.
data( ens.gene.ann.hg19 )
data( ens.gene.ann.hg19 )
A data.frame
Dataframe with 5 columns:
Gene: ensembl gene name
EnsID: ensembl gene id
Chromosome: Genomic Chromosome
Start: Genomic start of gene
End: Genomic end of gene
data.frame
A dataset containing the gene locations based on human genome reference hg38 that was obtained from BioMart.
data( ens.gene.ann.hg38 )
data( ens.gene.ann.hg38 )
A data.frame
Dataframe with 5 columns:
Gene: ensembl gene name
EnsID: ensembl gene id
Chromosome: Genomic Chromosome
Start: Genomic start of gene
End: Genomic end of gene
data.frame
Access Object featureChromosomes. This method returns a vector with feature chromosomes.
## S4 method for signature 'CopyNumberBreakPoints' featureChromosomes(object)
## S4 method for signature 'CopyNumberBreakPoints' featureChromosomes(object)
object |
An object of class |
a vector with feature chromosomes
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) featureChromosomes( bp )
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) featureChromosomes( bp )
Access Options featureInfo. This method returns a dataframe with feature annotations.
## S4 method for signature 'CopyNumberBreakPoints' featureInfo(object)
## S4 method for signature 'CopyNumberBreakPoints' featureInfo(object)
object |
of class |
data.frame with feature annotations
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) featureInfo( bp )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) featureInfo( bp )
Access Object featuresPerGene. This method returns a vector with gene-related features for a particular gene.
## S4 method for signature 'CopyNumberBreakPointGenes' featuresPerGene(object, geneName = NULL)
## S4 method for signature 'CopyNumberBreakPointGenes' featuresPerGene(object, geneName = NULL)
object |
An object of class |
geneName |
Exact Gene name as in the annotation |
a vector with gene-related features
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) featuresPerGene( bp, geneName="PCMTD2" )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) featuresPerGene( bp, geneName="PCMTD2" )
The GeneBreak package performs cohort based recurrent
gene breakpoint detection on copynumber data. It is possible
to use the output of the function CGHcall
from
the package CGHcall
or the function callBins
from the
package QDNAseq
as the input for this analysis.
Analysis starts with the function getBreakpoints
and continues with:bpFilter
to exclude certain breakpoints from the analysisaddGeneAnnotation
to add gene location informationbpGenes
to determine which features (probes/bins) are related to which genesbpStats
to determine which gene breaks are recurrent in the cohort
Access Object geneChromosomes. This method returns a vector with gene chromosomes.
## S4 method for signature 'CopyNumberBreakPointGenes' geneChromosomes(object)
## S4 method for signature 'CopyNumberBreakPointGenes' geneChromosomes(object)
object |
An object of class |
vector with gene chromosomes
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) geneChromosomes( bp )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) geneChromosomes( bp )
Access Options geneInfo. This method returns a dataframe with gene annotations.
## S4 method for signature 'CopyNumberBreakPointGenes' geneInfo(object)
## S4 method for signature 'CopyNumberBreakPointGenes' geneInfo(object)
object |
of class |
data.frame with gene annotations
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) geneInfo( bp )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) geneInfo( bp )
Builds the CopyNumberBreakPoints object from copynumber data and detects breakpoint locations.
getBreakpoints(data, data2 = NULL, first.rm = TRUE)
getBreakpoints(data, data2 = NULL, first.rm = TRUE)
data |
An object of class cghCall or an object of class QDNAseqCopyNumbers or a data.frame containing feature annotations ("Chromosome", "Start", "End", "FeatureName") followed by copy number segment values (rows are features, columns are subjects). |
data2 |
A "data.frame" containing copy number calls following feature annotations with the four columns ("Chromosome", "Start", "End", "FeatureName", ...). This is optional and allows CNA-associated breakpoint filtering. (see ?bpFilter) |
first.rm |
Remove the first 'artificial' breakpoint of the first DNA segment for each chromosome (default: first.rm=TRUE) |
The accuracy of chromosomal breakpoint locations depends on the quality and genomic resolution of processed copy number data. For CNA input data, we recommend to use established computational methods for CNA detection such as 'CGHcall' (Van De Wiel et al., 2007) for array-CGH or 'QDNAseq' (Scheinin et al., 2014) for MPS data, which both use the implemented Circular Binary Segmentation algorithm (Olshen et al. 2004).
Returns an object of class CopyNumberBreakPoints.
Van De Wiel, M.A. et al. (2007) CGHcall: calling aberrations for array CGH tumor profiles. Bioinformatics, 23, 892-894.
Scheinin, I. et al. (2014) DNA copy number analysis of fresh and formalin-fixed specimens by shallow whole-genome sequencing with identification and exclusion of problematic regions in the genome assembly. Genome Research, 24, 2022-2032.
Olshen, A.B. et al. (2004) Circular binary segmentation for the analysis of array-based DNA copy number data. 5, 557-572.
data( copynumber.data.chr20 ) breakpoints <- getBreakpoints( data = copynumber.data.chr20 ) ## or alternatively library(CGHcall) cgh <- copynumber.data.chr20 segmented <- data.frame( Chromosome=chromosomes(cgh), Start=bpstart(cgh), End=bpend(cgh), FeatureName=rownames(cgh), segmented(cgh)) called <- data.frame( Chromosome=chromosomes(cgh), Start=bpstart(cgh), End=bpend(cgh), FeatureName=rownames(cgh), calls(cgh)) breakpoints <- getBreakpoints( data = segmented, data2 = called ) ## options to inspect the data breakpoints accessOptions( breakpoints )
data( copynumber.data.chr20 ) breakpoints <- getBreakpoints( data = copynumber.data.chr20 ) ## or alternatively library(CGHcall) cgh <- copynumber.data.chr20 segmented <- data.frame( Chromosome=chromosomes(cgh), Start=bpstart(cgh), End=bpend(cgh), FeatureName=rownames(cgh), segmented(cgh)) called <- data.frame( Chromosome=chromosomes(cgh), Start=bpstart(cgh), End=bpend(cgh), FeatureName=rownames(cgh), calls(cgh)) breakpoints <- getBreakpoints( data = segmented, data2 = called ) ## options to inspect the data breakpoints accessOptions( breakpoints )
Access Object namesFeatures. This method returns a vector with feature names.
## S4 method for signature 'CopyNumberBreakPoints' namesFeatures(object)
## S4 method for signature 'CopyNumberBreakPoints' namesFeatures(object)
object |
An object of class |
a vector with feature names
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) namesFeatures( bp )
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) namesFeatures( bp )
Access Options recurrentGenes. This method returns a dataframe that contains genes that are recurrently affected across samples based on a FDR threshold.
## S4 method for signature 'CopyNumberBreakPointGenes' recurrentGenes(object, fdr.threshold = 0.1, summarize = TRUE, order.column = "FDR")
## S4 method for signature 'CopyNumberBreakPointGenes' recurrentGenes(object, fdr.threshold = 0.1, summarize = TRUE, order.column = "FDR")
object |
Output of bpStats(): an object of class |
fdr.threshold |
A numeric Genes with lower FDR are returned |
summarize |
A logical to determine whether to only output a selection of columns |
order.column |
Name of the column to sort output on |
data.frame with genes recurrently affected by breakpoints
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) bp <- bpStats( bp ) recurrentGenes( bp )
data( copynumber.data.chr20 ) data( ens.gene.ann.hg18 ) bp <- getBreakpoints( copynumber.data.chr20 ) bp <- bpFilter( bp ) bp <- addGeneAnnotation( bp, ens.gene.ann.hg18 ) bp <- bpGenes( bp ) bp <- bpStats( bp ) recurrentGenes( bp )
Access Object sampleNames. This method returns a vector with sample names.
## S4 method for signature 'CopyNumberBreakPoints' sampleNames(object)
## S4 method for signature 'CopyNumberBreakPoints' sampleNames(object)
object |
An object of class |
a vector with sample names
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) sampleNames( bp )
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) sampleNames( bp )
Access Object segmentData. This method returns a dataframe with segment values.
## S4 method for signature 'CopyNumberBreakPoints' segmentData(object)
## S4 method for signature 'CopyNumberBreakPoints' segmentData(object)
object |
An object of class |
a dataframe with segment values
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) segmentData( bp )
data( copynumber.data.chr20 ) bp <- getBreakpoints( copynumber.data.chr20 ) segmentData( bp )