Title: | A normalization method for Copy Number Aberration in cancer samples |
---|---|
Description: | Performs ratio, GC content correction and normalization of data obtained using low coverage (one read every 100-10,000 bp) high troughput sequencing. It performs a "discrete" normalization looking for the ploidy of the genome. It will also provide tumour content if at least two ploidy states can be found. |
Authors: | Stefano Berri <[email protected]>, Henry M. Wood <[email protected]>, Arief Gusnanto <[email protected]> |
Maintainer: | Stefano Berri <[email protected]> |
License: | GPL-2 |
Version: | 1.53.0 |
Built: | 2024-10-30 05:15:20 UTC |
Source: | https://github.com/bioc/CNAnorm |
addSmooth
segment ratio values in Package ‘CNAnorm’ using DNACopy
## S4 method for signature 'CNAnorm' addDNACopy(object, independent.arms = FALSE, ideograms = NULL, DNAcopy.smooth = list(), DNAcopy.weight = character(), DNAcopy.segment = list())
## S4 method for signature 'CNAnorm' addDNACopy(object, independent.arms = FALSE, ideograms = NULL, DNAcopy.smooth = list(), DNAcopy.weight = character(), DNAcopy.segment = list())
object |
An object of Class |
independent.arms |
Boolean. If TRUE chromosomes arms will be treated as independent, ideograms must be provided |
ideograms |
A data frame containing information about banding. See ?hg19_hs_ideogr for more information |
DNAcopy.smooth |
A list of parameters to be passed to function 'smooth.CNA' in package DNAcopy |
DNAcopy.weight |
A character value of one of these values. 'poisson', 'gaussian', 'twoTailQuantile', 'oneTailQuantile'. It specifies a way to give weight to each window depending on how much coverage in the normal deviates from the median for that chromosome. Options are listed in decreasing order of stringency. See Details |
DNAcopy.segment |
A list of parameters to be passed to function 'segment' in package DNAcopy. Parameters 'weights' and 'verbose' are not accepted |
'poisson': windows with coverage more or less than 2*sqrt(mean) from the mean are weighted down. The most stringent. Recommended for unbiased genome wide sequencing.
'gaussian': windows with coverage more or less than 2*sd from the median are weighted down. Recommended for genome wide sequencing where coverage in the normal is far from poisson distribution.
'twoTailQuantile': windows with coverage outside 5-95th quantile are weighted down. Recommended when coverage is far from a normal distribution - such as capture experiments -
'oneTailQuantile': windows with coverage lower than 5th quantile are weighted down. Recommended when coverage is far from a normal distribution - such as capture experiments. Does not weight down windows with very high coverage.
An object of class "CNAnorm"
signature(object = "CNAnorm")
Segment ratio values
on an object of class "CNAnorm"
. Returns the same object with
extra slots (segMean
, segID
)
Stefano Berri [email protected] and Arief Gusnanto [email protected]
Venkatraman, E. S. and Olshen, A. B. (2007) A faster circular binary segmentation algorithm for the analysis of array CGH data. Bioinformatics
segMean
, CNAnorm-class
, DNAcopy
, hg19_hs_ideogr
data(LS041) CN <- dataFrame2object(LS041) CN <- addDNACopy(CN)
data(LS041) CN <- dataFrame2object(LS041) CN <- addDNACopy(CN)
addSmooth
segment and smooth perform ratio values in Package ‘CNAnorm’
## S4 method for signature 'CNAnorm' addSmooth(object, lambda = 7, ...)
## S4 method for signature 'CNAnorm' addSmooth(object, lambda = 7, ...)
object |
An object of Class |
lambda |
Degree of smoothness. See reference for more details |
... |
Further arguments to pass to the function |
An object of class "CNAnorm"
signature(object = "CNAnorm")
Segment and smooth perform ratio values
on an object of class "CNAnorm"
. Returns the same object with
extra slot (ratio.s
)
Stefano Berri [email protected] and Arief Gusnanto [email protected]
Huang, J., Gusnanto, A., O'Sullivan, K., Staaf, J., Borg, A. and Pawitan, Y. (2007) Robust smooth segmentation approach for array CGH data analysis. Bioinformatics
data(LS041) CN <- dataFrame2object(LS041) CN.gcNorm <- gcNorm(CN, exclude = c("chrX", "chrY", "chrM")) CN.smooth <- addSmooth(CN)
data(LS041) CN <- dataFrame2object(LS041) CN.gcNorm <- gcNorm(CN, exclude = c("chrX", "chrY", "chrM")) CN.smooth <- addSmooth(CN)
chrs
returns/set the name of chromosomes/contigs
arms
retruns the name of the chromosome and arm. A data frame containing
ideogram information has to be provided. See ?hg19_hs_ideogr for an example
pos
returns/set the position of starting window. Be careful!
If you need to change data, it is better to change the input data and start over.
chrs(object) ## S4 method for signature 'CNAnorm' pos(object, show = "start") ## S4 method for signature 'CNAnorm' arms(object, banding_df)
chrs(object) ## S4 method for signature 'CNAnorm' pos(object, show = "start") ## S4 method for signature 'CNAnorm' arms(object, banding_df)
object |
An object of Class |
show |
The position to show: 'start', 'end' |
banding_df |
A data frame with infromation about ideogram |
chrs
and arms
return a character vector, pos
returns a numeric vector
Stefano Berri <[email protected]>
gcNorm
, CNAnorm-class
, hg19_hs_ideogr
data(LS041) data(hg19_hs_ideogr) CN <- dataFrame2object(LS041) dataFrameNames <- as.character(LS041$Chr) objectNames <- chrs(CN) armNames <- arms(CN, hg19_hs_ideogr) # check the names are, indeed, the same all(dataFrameNames == objectNames) # make shorter names, drop the first three letters ('chr') shortNames <- substr(chrs(CN),4,nchar(chrs(CN))) chrs(CN) <- shortNames # retrieve all new names unique(chrs(CN)) unique(armNames)
data(LS041) data(hg19_hs_ideogr) CN <- dataFrame2object(LS041) dataFrameNames <- as.character(LS041$Chr) objectNames <- chrs(CN) armNames <- arms(CN, hg19_hs_ideogr) # check the names are, indeed, the same all(dataFrameNames == objectNames) # make shorter names, drop the first three letters ('chr') shortNames <- substr(chrs(CN),4,nchar(chrs(CN))) chrs(CN) <- shortNames # retrieve all new names unique(chrs(CN)) unique(armNames)
This data is to provide an object to use in several examples without having to wait for computing it. To see how it was generated, see documentation of function peakPloidy.
data(CN)
data(CN)
A CNAnorm object
Class to Contain and Describe copy number aberration (CNA) data from low coverage (approx 0.01 - 0.5X) Next Generation Sequencing
Objects can be created by calls of the form new("CNAnorm", InData)
.
InData
:Object of class "InData"
. Contains input data
provided by the user. All slots have same length. Each element describe one window.
See Class "InData"
DerivData
:Object of class "DerivData"
. Contains data derived
from "InData"
. It can be Retrieved by the user, but methods should be used to populate
"DerivData"
. All slots have same length as input data. See Class DerivData
Res
:Object of class "Res"
. Contains slots with obtained results. See Class "Res"
Params
:Object of class "Params"
. Contains crucial parameters passed to some
of the methods for reusing in later steps or for documentation.
Summary of methods for class "CNAnorm"
. Type "methods ? methodName"
for more details about methodName.
signature(object = "CNAnorm")
: Retrieve Chromosomes/contig name
signature(object = "CNAnorm")
: Set Chromosomes/contig name
signature(object = "CNAnorm")
: Estimate ploidy
of the sample, tumor content and other results saved in Slot "Res"
signature(x = "CNAnorm")
: Returns number of element/windows
signature(x = "CNAnorm")
: Produce on object of class "CNAnorm"
with a subser of windows
signature(object = "CNAnorm")
: Plot annotated normalized genome copy number
signature(object = "CNAnorm")
: Plot peaks and estimated/validated ploidy
signature(object = "CNAnorm")
: Retrieve Chromosomes/contig position
signature(object = "CNAnorm")
: Set Chromosomes/contig position
signature(object = "CNAnorm")
: Retrieve ratio (Test/Control).
If gcNorm was called, ratio is GC normalized
signature(object = "CNAnorm")
: Retrieve normalized ratio (not smoothed)
signature(object = "CNAnorm")
: Retrieve smoothed ratio
signature(object = "CNAnorm")
: Retrieve normalized smoothed ratio
signature(object = "CNAnorm")
: Retrieve segmented ratio (as provided by DNAcopy)
signature(object = "CNAnorm")
: Retrieve normalized segmented ratio
Stefano Berri [email protected] and Arief Gusnanto [email protected]
CNA-norm: Discrete Normalization of Copy Number Alteration data from clinical samples (in preparation)
InData
, DerivData
for documentation on the slots.
data(LS041) CNA <- new("CNAnorm", InData = new("InData", Chr = as.character(LS041$Chr), Pos = LS041$Pos, Test = LS041$Test, Norm = LS041$Norm, GC = LS041$GC))
data(LS041) CNA <- new("CNAnorm", InData = new("InData", Chr = as.character(LS041$Chr), Pos = LS041$Pos, Test = LS041$Test, Norm = LS041$Norm, GC = LS041$GC))
"CNAnorm"
workflowThis function is a wrapper to use for a fully automated CNAnorm workflow where interactivity is not required. It contains MOST possible paramenters. Defaults are set to to run a standard and conservative workflow.
CNAnormWorkflow(dataFrame, gc.do=FALSE, gc.exclude=character(0), gc.maxNumPoints=10000, smooth.do=TRUE, smooth.lambda=7, smooth.other=list(), peak.method='closest', peak.exclude=character(0), peak.ploidyToTest=12, peak.sd=5, peak.dThresh=0.01, peak.n=2048, peak.adjust=.9, peak.force.smooth=TRUE, peak.reg=FALSE, peak.ds=1.5, peak.zero.count=FALSE, peak.other=list(), DNAcopy.do=TRUE, DNAcopy.independent.arms=FALSE, DNAcopy.ideograms=NULL, DNAcopy.smooth=list(), DNAcopy.segment=list(), DNAcopy.weight=character(), dNorm.normBy=NULL)
CNAnormWorkflow(dataFrame, gc.do=FALSE, gc.exclude=character(0), gc.maxNumPoints=10000, smooth.do=TRUE, smooth.lambda=7, smooth.other=list(), peak.method='closest', peak.exclude=character(0), peak.ploidyToTest=12, peak.sd=5, peak.dThresh=0.01, peak.n=2048, peak.adjust=.9, peak.force.smooth=TRUE, peak.reg=FALSE, peak.ds=1.5, peak.zero.count=FALSE, peak.other=list(), DNAcopy.do=TRUE, DNAcopy.independent.arms=FALSE, DNAcopy.ideograms=NULL, DNAcopy.smooth=list(), DNAcopy.segment=list(), DNAcopy.weight=character(), dNorm.normBy=NULL)
All arguments are explained in the relative functions
dataFrame |
A data frame with columns Chr, Pos, Test, Norm and optional GC. See dataFrame2object |
gc.do |
Specify if GC correction need to be done. See gcNorm |
gc.exclude |
See gcNorm |
gc.maxNumPoints |
See gcNorm |
smooth.do |
Specify if smoothing need to be done. See addSmooth |
smooth.lambda |
See addSmooth |
smooth.other |
A list of other parameters to pass to the smoothing function. See addSmooth |
peak.method |
See peakPloidy |
peak.exclude |
See peakPloidy |
peak.ploidyToTest |
See peakPloidy |
peak.sd |
See peakPloidy |
peak.dThresh |
See peakPloidy |
peak.n |
See peakPloidy |
peak.adjust |
See peakPloidy |
peak.force.smooth |
See peakPloidy |
peak.reg |
See peakPloidy |
peak.ds |
See peakPloidy |
peak.zero.count |
See peakPloidy |
peak.other |
A list of other parameters to be passed to funtions for peak detection See peakPloidy |
DNAcopy.do |
Specify if segmentation with DNAcopy need to be done. See addDNACopy |
DNAcopy.independent.arms |
See addDNACopy |
DNAcopy.ideograms |
See addDNACopy |
DNAcopy.smooth |
See addDNACopy |
DNAcopy.segment |
See addDNACopy |
DNAcopy.weight |
See addDNACopy |
dNorm.normBy |
See discreteNorm |
An object of Class "CNAnorm"
Stefano Berri <[email protected]>
dataFrame2object
, gcNorm
, addSmooth
,
peakPloidy
, addDNACopy
, discreteNorm
data(LS041) CN <- CNAnormWorkflow(LS041)
data(LS041) CN <- CNAnormWorkflow(LS041)
"CNAnorm"
Convert a data frame with column: Chr, Pos, Test, Norm and
optional GC into object of class "CNAnorm"
dataFrame2object(dataFrame)
dataFrame2object(dataFrame)
dataFrame |
A data frame with columns Chr, Pos, Test, Norm and optional GC |
An object of Class "CNAnorm"
Stefano Berri <[email protected]>
CNAnorm-class
, InData-class
, data.frame
data(LS041) CN <- dataFrame2object(LS041)
data(LS041) CN <- dataFrame2object(LS041)
A Class containing data derived from InData used for further computation and plotting.
Objects can be created by calls of the form new("DerivData")
, however DerivData is typically
populated using methods. If a slot has not been populated yet, it has zero length, otherwise slots
have the same length as InData
.
ratio
:Numeric vector with ratio Test/Normal. Optionally GC corrected.
ratio.s
:Numeric vector with smoothed ratio.
ratio.n
:Numeric vector with normalized ratio.
ratio.s.n
:Numeric vector with normalized and smoothed ratio.
segID
:Numeric vecotr with ID of segmented data (as provided by DNACopy).
Each segment
has a different ID.
segMean
:Numeric vector with mean value of the segment (as provided by DNACopy.)
segMean.n
:Numeric vector with normalized segMean.
ok4density
:Logical vector. Specify wich values have been used to calculate density.
signature(x = "DerivData")
: Returns number of windows.
Stefano Berri and Arief Gusnanto
Gusnanto, A., Wood, H.M., Pawitan, Y., Rabbitts, P. and Berri, S. (2011) Correcting for cancer genome size and tumor cell content enables better estimation of copy number alterations from next generation sequence data. Bioinformatics
data(LS041) inObject <- new("InData", Chr = as.character(LS041$Chr), Pos = LS041$Pos, Test = LS041$Test, Norm = LS041$Norm, GC = LS041$GC) CNA <- new("CNAnorm", InData = inObject)
data(LS041) inObject <- new("InData", Chr = as.character(LS041$Chr), Pos = LS041$Pos, Test = LS041$Test, Norm = LS041$Norm, GC = LS041$GC) CNA <- new("CNAnorm", InData = inObject)
discreteNorm
performs normalization of data using information on ploidy.
Implicitly calls "validation"
if no validation was performed
## S4 method for signature 'CNAnorm' discreteNorm(object, normBy = object)
## S4 method for signature 'CNAnorm' discreteNorm(object, normBy = object)
object |
An object of Class |
normBy |
An object of Class |
An object of class "CNAnorm"
Stefano Berri [email protected] and Arief Gusnanto [email protected]
Gusnanto, A., Wood, H.M., Pawitan, Y., Rabbitts, P. and Berri, S. (2011) Correcting for cancer genome size and tumour cell content enables better estimation of copy number alterations from next generation sequence data. Bioinformatics
data(CN) # see peakPloidy documentation to know how object CN is created CN <- discreteNorm(CN)
data(CN) # see peakPloidy documentation to know how object CN is created CN <- discreteNorm(CN)
exportTable
write a table with normalised values of each window. A wrapper to "write.table"
## S4 method for signature 'CNAnorm' exportTable(object, file = "CNAnorm_table.tab", show = 'ratio', sep = "\t", row.names = FALSE, ...)
## S4 method for signature 'CNAnorm' exportTable(object, file = "CNAnorm_table.tab", show = 'ratio', sep = "\t", row.names = FALSE, ...)
object |
an object of Class |
file |
name of the file to save to |
show |
what should be reported in the table: |
sep |
the field separator string. |
row.names |
either a logical value indicating whether the row number should be written or a character vector of row names to be written. |
... |
Extra arguments to be passed to |
It produces a tab delimited text file with the following columns:
Chr: Chromosome/contig name.
Pos: Starting position of the window.
Ratio: Ratio Test/Normal for each window after GC correction.
Ratio.n: Ratio Test/Normal or ploidy for each window after normalisation.
Ratio.s.n: Smoothed and normalised ratio Test/Normal or ploidy for each window.
SegMean: Mean of the segment this window belongs to.
SegMean.n: Normalised mean ratio Test/Normal or ploidy of the segment this window belongs to.
An object of class "CNAnorm"
Stefano Berri [email protected]
data(CN) CN <- validation(CN) CN <- discreteNorm(CN) exportTable(CN, file = "CNAnorm_table.tab", show = 'ploidy')
data(CN) CN <- validation(CN) CN <- discreteNorm(CN) exportTable(CN, file = "CNAnorm_table.tab", show = 'ploidy')
gcNorm
perform GC content normalization on ratio Test/Normal in Package ‘CNAnorm’
## S4 method for signature 'CNAnorm' gcNorm(object, exclude = character(0), maxNumPoints = 10000)
## S4 method for signature 'CNAnorm' gcNorm(object, exclude = character(0), maxNumPoints = 10000)
object |
An object of Class |
exclude |
A character vector with name of chromosomes/contigues not to use to calculate GC content correction. All genome, however, will be corrected |
maxNumPoints |
Maximum number of data points to fit the loess correction.
For computational pourposes, if the number of points in |
An object of class "CNAnorm"
signature(object = "CNAnorm")
Perform GC content correction on
an object of class "CNAnorm"
. Returns the same object with corrected ratio
Stefano Berri <[email protected]>
data(LS041) CN <- dataFrame2object(LS041) # correct for GC content, but ignoring data from sex chromosomes and # mitocondria CN.gcNorm <- gcNorm(CN, exclude = c("chrX", "chrY", "chrM"))
data(LS041) CN <- dataFrame2object(LS041) # correct for GC content, but ignoring data from sex chromosomes and # mitocondria CN.gcNorm <- gcNorm(CN, exclude = c("chrX", "chrY", "chrM"))
This data object is used by some plotting methods and contains the default values. User can change graphical parameters by changing this object
The object consists of several layers refering to different plots and different properties. Here an indicative description
gPar$genome: prameters here refer to the plot produced by plotGenome
graphical parameters: see ?par
$colors: specify colors $cex: specify relative size - points, text... $lwd: specify line width $lty: specify line type - solid, dashed $mar: specify margins
data(gPar)
data(gPar)
A S3 object
This is bundles data that can be provided to plotGenome in order to plot location of the centromere. In future release it might be used to produce an ideogram too
data(hg19_hs_ideogr)
data(hg19_hs_ideogr)
A data.frame
A Class containing input data for CNA
Objects can be created by calls of the form new("InData", Chr, Pos, Test, Norm, ratio, GC)
.
Chr
:Object of class "character"
. Name of the Chromosomes/Contigs of each window.
Pos
:Object of class "numeric"
. Starting position of the each window.
Test
:Object of class "numeric"
. Number of reads from Test in each window.
Norm
:Object of class "numeric"
. Number of reads from Normal (Control) in each window.
ratio
:Object of class "numeric"
. Ratio Test/Control
in each window. Automatically computed if Test and Norm are provided, or
user generated if Test and Norm are not know.
GC
:Object of class "numeric"
. GC content of each window.
signature(x = "InData")
: Returns number of windows from input data.
Stefano Berri
Gusnanto, A., Wood, H.M., Pawitan, Y., Rabbitts, P. and Berri, S. (2011) Correcting for cancer genome size and tumor cell content enables better estimation of copy number alterations from next generation sequence data. Bioinformatics
data(LS041) inObject <- new("InData", Chr = as.character(LS041$Chr), Pos = LS041$Pos, Test = LS041$Test, Norm = LS041$Norm, GC = LS041$GC) CNA <- new("CNAnorm", InData = inObject)
data(LS041) inObject <- new("InData", Chr = as.character(LS041$Chr), Pos = LS041$Pos, Test = LS041$Test, Norm = LS041$Norm, GC = LS041$GC) CNA <- new("CNAnorm", InData = inObject)
This data set provide reads in tumor and matched blood for patient LS041. Each row has information about non-overlapping window across the genome. In particular it reports: chromosome name, starting position of the window (1 based), number of mapped reads in the test (lung tumor), number of reads in the control (matched blood) and GC content of the window.
data(LS041)
data(LS041)
A dataframe
Gusnanto, A., Wood, H.M., Pawitan, Y., Rabbitts, P. and Berri, S. (2011) Correcting for cancer genome size and tumour cell content enables better estimation of copy number alterations from next generation sequence data. Bioinformatics
A Class containing some Parameters used in the analysis. Not heavily used at the moment.
Objects can be created by calls of the form new("Params")
, it is usually
iniziated and populated with methods and functions of class CNAnorm
.
method
:variable of class "character"
. Record if the peakPloidy
function was called using density
or mixture
.
density.n
:The variable "n"
used when calling peakPloidy
.
This variable is saved so that can be used later for drawing plots.
density.adjust
:The variable "adjust"
used when calling peakPloidy
.
This variable is saved so that can be used later for drawing plots
gc.excludeFromGCNorm
:Vector of class "character"
.
Name of the Chromosomes/Contigs not used for GC content correction.
gc.maxNumPoints
:One element vector of class "numeric"
.
Specify how many points to use for GC correction
gp.excludeFromDensity
:Vector of class "character"
.
Name of the Chromosomes/Contigs not used for peak guessing
signature(x = "Params")
Stefano Berri
Gusnanto, A., Wood, H.M., Pawitan, Y., Rabbitts, P. and Berri, S. (2011) Correcting for cancer genome size and tumor cell content enables better estimation of copy number alterations from next generation sequence data. Bioinformatics
data(LS041) inObject <- new("InData", Chr = as.character(LS041$Chr), Pos = LS041$Pos, Test = LS041$Test, Norm = LS041$Norm, GC = LS041$GC) CNA <- new("CNAnorm", InData = inObject)
data(LS041) inObject <- new("InData", Chr = as.character(LS041$Chr), Pos = LS041$Pos, Test = LS041$Test, Norm = LS041$Norm, GC = LS041$GC) CNA <- new("CNAnorm", InData = inObject)
peakPloidy
Estimate most likely ploidy of genome looking at distribution of smoothed ratio.
## S4 method for signature 'CNAnorm' peakPloidy(object, method = 'mixture', exclude = character(0), ploidyToTest = 12, sd = 5, dThresh = 0.01, n = 2048, adjust = .9, force.smooth = TRUE, reg = FALSE, ds = 1.5, zero.cont = FALSE, ...)
## S4 method for signature 'CNAnorm' peakPloidy(object, method = 'mixture', exclude = character(0), ploidyToTest = 12, sd = 5, dThresh = 0.01, n = 2048, adjust = .9, force.smooth = TRUE, reg = FALSE, ds = 1.5, zero.cont = FALSE, ...)
object |
An object of Class |
exclude |
A character vector with names of Chromosomes/Contigs not to use to estimate ploidy. |
method |
A character element matching either |
ploidyToTest |
Maximum ploidy allowed. Warnings! Computation time
increases exponentially with this parameter if using |
adjust |
The |
n |
The |
force.smooth |
If the input object does not have smoothed ratio, it will smooth
using |
sd |
Parameter to filter outliers. Values greater than i median + sd * standard deivations will be ignored while detecting peaks and ploidy. |
dThresh |
Parameter to filter outliers. Values with a density lower than max(density)*dThresh will be ignored while detecting peaks and ploidy. |
reg |
Parameter for mixture model: If set TRUE, the starting point for EM algorithm will be optimized through several methods including regular grid on the ratio distribution. The default is FALSE, by which the starting values are taken from the quantiles of the distribution. |
ds |
Parameter for mixture model: A constraint in EM algorithm of minimum distance between mean estimates, in terms of median standard deviation of the mixture components. |
zero.cont |
Parameter for mixture model: An argument for the mixture model. If set TRUE, the EM algorithm considers exactly-zero ratios as a mixture component. |
... |
Extra parameters to be passed to funtions for peak detection, specific to each of the methods (deinsity or mixture), se below for details. |
An object of class "CNAnorm"
Other optional parameters to be passed (...)
mixture method
density method
peakRatioThreshold used to call a peak. Peaks smaller than
maxPeakHight/peakRatio
are not considered peaks.
spacingToleranceA parameter smaller than 1 which describes how strict the program should be on alternative solutions. Only solution with the best R^2 >= max(R^2)*spacingTolerance will be considered as valid.
interceptRatioMinimum value for the intercept of the linear model. Ideally, should be zero, but the default allows a little flexibility.
Stefano Berri [email protected] and Arief Gusnanto [email protected]
Gusnanto, A., Wood, H.M., Pawitan, Y., Rabbitts, P. and Berri, S. (2011) Correcting for cancer genome size and tumour cell content enables better estimation of copy number alterations from next generation sequence data. Bioinformatics
data(LS041) CN <- dataFrame2object(LS041) chr2skip <- c("chrY", "chrM") CN <- gcNorm(CN, exclude = chr2skip) CN <- addSmooth(CN, lambda = 3) CN <- peakPloidy(CN, exclude = chr2skip) # this object CN is what you obtain when you load # data(CN)
data(LS041) CN <- dataFrame2object(LS041) chr2skip <- c("chrY", "chrM") CN <- gcNorm(CN, exclude = chr2skip) CN <- addSmooth(CN, lambda = 3) CN <- peakPloidy(CN, exclude = chr2skip) # this object CN is what you obtain when you load # data(CN)
plotGenome
plot normalized ratio and optionally segmented and/or smoothed
normalized ratio values in Package ‘CNAnorm’. It also shows annotation.
## S4 method for signature 'CNAnorm' plotGenome(object, maxRatio = 8, minRatio = -1, superimpose = character(0), gPar = NULL, numHorLables = 10, colorful = FALSE, show.centromeres = TRUE, idiogram = NULL, fixVAxes = FALSE, supLineColor = character(0), supLineCex = character(0), dot.cex = .2, ...)
## S4 method for signature 'CNAnorm' plotGenome(object, maxRatio = 8, minRatio = -1, superimpose = character(0), gPar = NULL, numHorLables = 10, colorful = FALSE, show.centromeres = TRUE, idiogram = NULL, fixVAxes = FALSE, supLineColor = character(0), supLineCex = character(0), dot.cex = .2, ...)
object |
An object of Class |
maxRatio |
The maximum ratio to be shown on the plot. Values or ratio greater than maxRatio will be displayed as green triangulars |
minRatio |
The minimum ratio to be shown on the plot. Values or ratio smaller than minRatio will be displayed as green triangulars |
superimpose |
A character verctor with one or both of the following:
|
numHorLables |
. Number of maximum horizontal lables. The function will try to annotate numHorLables so that they are approximately equally spaced. |
colorful |
A switch to decide if the background dots representing the ratio of each window should be gray or colored according their value in relation to the peak closest to the median |
show.centromeres |
A switch to decide if location of centromere are displayed on the graph. The location of the centromere is stored in idiogram |
idiogram |
A data frame containing banding information. if NULL -default- human information will be loaded by data(hg19_hs_ideogr) |
fixVAxes |
A switch to decide if the vertical axes should be fixed to minRatio and maxRatio or fit the data within minRatio and maxRatio. |
gPar |
a S3 object with all graphical parameters. If NULL (default) data(gPar) is called |
supLineColor |
A three element character vector with colors to be used for
|
supLineCex |
A two element vector with |
dot.cex |
size of the dots in the plot |
... |
Further arguments to pass to the function |
An object of class "CNAnorm"
Stefano Berri [email protected] and Arief Gusnanto [email protected]
plot
, par
, peakPloidy
,
gPar
, hg19_hs_ideogr
data(CN) # see peakPloidy documentation to know how object CN is created CN <- addDNACopy(CN) CN <- validation(CN) CN <- discreteNorm(CN) plotGenome(CN, superimpose = 'DNACopy')
data(CN) # see peakPloidy documentation to know how object CN is created CN <- addDNACopy(CN) CN <- validation(CN) CN <- discreteNorm(CN) plotGenome(CN, superimpose = 'DNACopy')
plotPeaks
plot annotated distribution of ratio Test/Normal
## S4 method for signature 'CNAnorm' plotPeaks(object, special1 = character(0), special2 = character(0), show ='suggested', ...)
## S4 method for signature 'CNAnorm' plotPeaks(object, special1 = character(0), special2 = character(0), show ='suggested', ...)
object |
An object of Class |
special1 |
The chromosome/contig whose distribution will be shown with a different color |
special2 |
The chromosome/contig whose distribution will be shown with a different color |
show |
A character verctor with one or both of the following:
|
... |
Further arguments to pass to the function |
Stefano Berri [email protected]
data(CN) # see peakPloidy documentation to know how object CN is created plotPeaks(CN, special1 = 'chrX', special2 = 'chrY')
data(CN) # see peakPloidy documentation to know how object CN is created plotPeaks(CN, special1 = 'chrX', special2 = 'chrY')
ratio
returns the Test/Normal ratio from an object of class CNAnorm
.
ratio is corrected for GC content if gcNorm
was called.
ratio.n
returns the Test/Normal normalized ratio from an object of class CNAnorm
after normalization. Its input is ratio(object)
ratio.s
returns the Test/Normal smoothed ratio from an object of class CNAnorm
Its input is ratio(object)
ratio.s.n
returns the Test/Normal smoothed and normalized ratio
from an object of class CNAnorm
. Its input is ratio.s(object)
segMean
returns the mean of the segments as produced by DNACopy
segMean.n
returns the normalized mean of the segments
ratio(object) ratio.n(object) ratio.s(object) ratio.s.n(object) segMean(object) segMean.n(object)
ratio(object) ratio.n(object) ratio.s(object) ratio.s.n(object) segMean(object) segMean.n(object)
object |
An object of Class |
A numeric vector
Stefano Berri <[email protected]>
gcNorm
, CNAnorm-class
, DNAcopy
data(LS041) CN <- dataFrame2object(LS041) ratio.original <- ratio(CN) CN.gcNorm <- gcNorm(CN, exclude = c("chrX", "chrY", "chrM")) ratio.gc.corrected <- ratio(CN.gcNorm)
data(LS041) CN <- dataFrame2object(LS041) ratio.original <- ratio(CN) CN.gcNorm <- gcNorm(CN, exclude = c("chrX", "chrY", "chrM")) ratio.gc.corrected <- ratio(CN.gcNorm)
sugg.peaks
returns the location of peaks before normalization as
suggested by peakPloidy
.
sugg.ploidy
returns the ploidy of the peaks as suggested by
peakPloidy
.
valid.peaks
returns the location of peaks before normalization as
validated after calling method "validation"
valid.ploidy
returns the validated ploidy of the peaks as validated
after calling method "validation"
sugg.peaks(object) sugg.ploidy(object) valid.peaks(object) valid.ploidy(object)
sugg.peaks(object) sugg.ploidy(object) valid.peaks(object) valid.ploidy(object)
object |
An object of Class |
A numeric vector
Stefano Berri <[email protected]>
gcNorm
, CNAnorm-class
, DNAcopy
data(CN) # see peakPloidy documentation to know how object CN is created plot(sugg.ploidy(CN), sugg.peaks(CN))
data(CN) # see peakPloidy documentation to know how object CN is created plot(sugg.ploidy(CN), sugg.peaks(CN))
validation
segment and smooth perform ratio values in Package ‘CNAnorm’
## S4 method for signature 'CNAnorm' validation(object, peaks = sugg.peaks(object), ploidy = sugg.ploidy(object))
## S4 method for signature 'CNAnorm' validation(object, peaks = sugg.peaks(object), ploidy = sugg.ploidy(object))
object |
An object of Class |
peaks |
The user validated location (ratio Test/Normal) of the peaks before normalization. |
ploidy |
The user validated ploidy of the peaks before normalization. |
An object of class "CNAnorm"
It is implicitly called by discreteNorm
if no validation was manually performed
Stefano Berri [email protected]
data(CN) # see peakPloidy documentation to know how object CN is created CN <- validation(CN)
data(CN) # see peakPloidy documentation to know how object CN is created CN <- validation(CN)