
Last chance! 50% off unlimited learning
Sale ends in
Segment total copy numbers and allele B fractions using the Paired PSCBS method [1]. This method requires matched normals. This is a low-level segmentation method. It is intended to be applied to one tumor-normal sample at the time.
# S3 method for default
segmentByPairedPSCBS(CT, thetaT=NULL, thetaN=NULL, betaT=NULL, betaN=NULL, muN=NULL,
rho=NULL, chromosome=0, x=NULL, alphaTCN=0.009, alphaDH=0.001, undoTCN=0, undoDH=0,
..., avgTCN=c("mean", "median"), avgDH=c("mean", "median"),
flavor=c("tcn&dh", "tcn,dh", "sqrt(tcn),dh", "sqrt(tcn)&dh", "tcn"), tbn=is.null(rho),
joinSegments=TRUE, knownSegments=NULL, dropMissingCT=TRUE, seed=NULL, verbose=FALSE,
preserveScale=FALSE)
Returns the segmentation results as a PairedPSCBS
object.
A numeric
vector
of J tumor total copy number (TCN)
ratios in [0,+Inf
) (due to noise, small negative values are
also allowed). The TCN ratios are typically scaled such that
copy-neutral diploid loci have a mean of two.
(alternative) As an alternative to specifying
tumor TCN ratios relative to the match normal by
argument CT
, on may specify total tumor and normal
signals seperately, in which case the TCN ratios CT
are
calculated as
A numeric
vector
of J tumor allele B fractions (BAFs)
in [0,1] (due to noise, values may be slightly outside as well)
or NA
for non-polymorphic loci.
A numeric
vector
of J matched normal BAFs in [0,1]
(due to noise, values may be slightly outside as well) or NA
for non-polymorphic loci.
An optional numeric
vector
of J genotype calls in
{0,1/2,1} for AA, AB, and BB, respectively,
and NA
for non-polymorphic loci.
If not given, they are estimated from the normal BAFs using
callNaiveGenotypes
as described in [2].
(alternative to betaT
and betaN
/muN
)
A numeric
vector
of J decrease-of-heterozygosity signals (DHs)
in [0,1] (due to noise, values may be slightly larger than one
as well). By definition, DH should be NA
for homozygous loci
and for non-polymorphic loci.
(Optional) An integer
scalar (or a vector
of length J),
which can be used to specify which chromosome each locus belongs to
in case multiple chromosomes are segments.
This argument is also used for annotation purposes.
Optional numeric
vector
of J genomic locations.
If NULL
, index locations 1:J
are used.
The significance levels for segmenting total copy numbers (TCNs) and decrease-in-heterozygosity signals (DHs), respectively.
Non-negative numeric
s. If greater than 0,
then a cleanup of segmentions post segmentation is done.
See argument undo
of segmentByCBS
() for more
details.
A character
string specifying how to calculating
segment mean levels after change points have been
identified.
Additional arguments passed to segmentByCBS
().
A character
specifying what type of segmentation and
calling algorithm to be used.
If TRUE
, betaT
is normalized before segmentation
using the TumorBoost method [2], otherwise not.
If TRUE
, there are no gaps between neighboring
segments.
If FALSE
, the boundaries of a segment are defined by the support
that the loci in the segments provides, i.e. there exist a locus
at each end point of each segment. This also means that there
is a gap between any neighboring segments, unless the change point
is in the middle of multiple loci with the same position.
The latter is what DNAcopy::segment()
returns.
Optional data.frame
specifying
non-overlapping known segments. These segments must
not share loci. See findLargeGaps
() and gapsToSegments
().
If TRUE
, loci for which 'CT' is missing
are dropped, otherwise not.
An (optional) integer
specifying the random seed to be
set before calling the segmentation method. The random seed is
set to its original state when exiting. If NULL
, it is not set.
See Verbose
.
Defunct - gives an error is specified.
The "DNAcopy::segment" implementation of CBS uses approximation through random sampling for some estimates. Because of this, repeated calls using the same signals may result in slightly different results, unless the random seed is set/fixed.
Although it is possible to segment each chromosome independently using Paired PSCBS, we strongly recommend to segment whole-genome (TCN,BAF) data at once. The reason for this is that downstream CN-state calling methods, such as the AB and the LOH callers, performs much better on whole-genome data. In fact, they may fail to provide valid calls if done chromosome by chromosome.
The total copy number signals as well as any optional positions
must not contain missing values, i.e. NA
s or NaN
s.
If there are any, an informative error is thrown.
Allele B fractions may contain missing values, because such are
interpreted as representing non-polymorphic loci.
None of the input signals may have infinite values, i.e. -Inf
or +Inf
.
If so, an informative error is thrown.
If allele B fractions for the matched normal (betaN
) are
not available, but genotypes (muN
) are, then it is possible
to run a version of Paired PSCBS where TumorBoost normalization
of the tumor allele B fractions is skipped. In order for this
to work, argument tbn
must be set to FALSE
.
Henrik Bengtsson
Internally segmentByCBS
() is used for segmentation.
The Paired PSCBS segmentation method does not support weights.
[1] A.B. Olshen, H. Bengtsson, P. Neuvial, P.T. Spellman, R.A. Olshen, V.E. Seshan, Parent-specific copy number in paired tumor-normal studies using circular binary segmentation, Bioinformatics, 2011
[2] H. Bengtsson, P. Neuvial and T.P. Speed, TumorBoost: Normalization of allele-specific tumor copy numbers from a single pair of tumor-normal genotyping microarrays, BMC Bioinformatics, 2010
Internally, callNaiveGenotypes
is used to
call naive genotypes, normalizeTumorBoost
is
used for TumorBoost normalization, and segmentByCBS
() is used
to segment TCN and DH separately.
To segment tumor total copy numbers and allele B fractions
without a matched normal, see segmentByNonPairedPSCBS
().
To segment total copy-numbers, or any other unimodal signals,
see segmentByCBS
().
verbose <- R.utils::Arguments$getVerbose(-10*interactive(), timestamp=TRUE)
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
# Load SNP microarray data
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
data <- PSCBS::exampleData("paired.chr01")
str(data)
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
# Paired PSCBS segmentation
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
# Drop single-locus outliers
dataS <- dropSegmentationOutliers(data)
# Speed up example by segmenting fewer loci
dataS <- dataS[seq(from=1, to=nrow(data), by=10),]
str(dataS)
R.oo::attachLocally(dataS)
# Paired PSCBS segmentation
fit <- segmentByPairedPSCBS(CT, betaT=betaT, betaN=betaN,
chromosome=chromosome, x=x,
seed=0xBEEF, verbose=verbose)
print(fit)
# Plot results
plotTracks(fit)
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
# Bootstrap segment level estimates
# (used by the AB caller, which, if skipped here,
# will do it automatically)
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
fit <- bootstrapTCNandDHByRegion(fit, B=100, verbose=verbose)
print(fit)
plotTracks(fit)
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
# Calling segments in allelic balance (AB)
# NOTE: Ideally, this should be done on whole-genome data
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
# Explicitly estimate the threshold in DH for calling AB
# (which be done by default by the caller, if skipped here)
deltaAB <- estimateDeltaAB(fit, flavor="qq(DH)", verbose=verbose)
print(deltaAB)
## [1] 0.1657131
fit <- callAB(fit, delta=deltaAB, verbose=verbose)
print(fit)
plotTracks(fit)
# Even if not explicitly specified, the estimated
# threshold parameter is returned by the caller
stopifnot(fit$params$deltaAB == deltaAB)
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
# Calling segments in loss-of-heterozygosity (LOH)
# NOTE: Ideally, this should be done on whole-genome data
# - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
# Explicitly estimate the threshold in C1 for calling LOH
# (which be done by default by the caller, if skipped here)
deltaLOH <- estimateDeltaLOH(fit, flavor="minC1|nonAB", verbose=verbose)
print(deltaLOH)
## [1] 0.625175
fit <- callLOH(fit, delta=deltaLOH, verbose=verbose)
print(fit)
plotTracks(fit)
# Even if not explicitly specified, the estimated
# threshold parameter is returned by the caller
stopifnot(fit$params$deltaLOH == deltaLOH)
Run the code above in your browser using DataLab