segmentationPSCBS: PSCBS segmentation

View source: R/segmentationPSCBS.R

segmentationPSCBSR Documentation

PSCBS segmentation

Description

Alternative segmentation function using the PSCBS package. This function is called via the fun.segmentation argument of runAbsoluteCN. The arguments are passed via args.segmentation.

Usage

segmentationPSCBS(
  normal,
  tumor,
  log.ratio,
  seg,
  plot.cnv,
  sampleid,
  weight.flag.pvalue = 0.01,
  alpha = 0.005,
  undo.SD = NULL,
  flavor = "tcn&dh",
  tauA = 0.03,
  vcf = NULL,
  tumor.id.in.vcf = 1,
  normal.id.in.vcf = NULL,
  max.segments = NULL,
  boost.on.target.max.size = 30,
  min.logr.sdev = 0.15,
  prune.hclust.h = NULL,
  prune.hclust.method = "ward.D",
  chr.hash = NULL,
  additional.cmd.args = "",
  centromeres = NULL,
  ...
)

Arguments

normal

Coverage data for normal sample. Ignored in this function.

tumor

Coverage data for tumor sample.

log.ratio

Copy number log-ratios, one for each exon in coverage file.

seg

If segmentation was provided by the user, this data structure will contain this segmentation. Useful for minimal segmentation functions. Otherwise PureCN will re-segment the data. This segmentation function ignores this user provided segmentation.

plot.cnv

Segmentation plots.

sampleid

Sample id, used in output files.

weight.flag.pvalue

Flag values with one-sided p-value smaller than this cutoff.

alpha

Alpha value for CBS, see documentation for the segment function.

undo.SD

undo.SD for CBS, see documentation of the segment function. If NULL, try to find a sensible default.

flavor

Flavor value for PSCBS. See segmentByNonPairedPSCBS.

tauA

tauA argument for PSCBS. See segmentByNonPairedPSCBS.

vcf

Optional VCF object with germline allelic ratios.

tumor.id.in.vcf

Id of tumor in case multiple samples are stored in VCF.

normal.id.in.vcf

Id of normal in in VCF. If NULL, use unpaired PSCBS.

max.segments

If not NULL, try a higher undo.SD parameter if number of segments exceeds the threshold.

boost.on.target.max.size

When off-target regions are noisy compared to on-target, try to find small segments of specified maximum size that might be missed to due the increased noise. Set to 0 to turn boosting off.

min.logr.sdev

Minimum log-ratio standard deviation used in the model. Useful to make fitting more robust to outliers in very clean data.

prune.hclust.h

Height in the hclust pruning step. Increasing this value will merge segments more aggressively. If NULL, try to find a sensible default.

prune.hclust.method

Cluster method used in the hclust pruning step. See documentation for the hclust function.

chr.hash

Mapping of non-numerical chromsome names to numerical names (e.g. chr1 to 1, chr2 to 2, etc.). If NULL, assume chromsomes are properly ordered.

additional.cmd.args

character(1). Ignored.

centromeres

A GRanges with centromere positions. If not NULL, add breakpoints at centromeres.

...

Additional parameters passed to the segmentByNonPairedPSCBS function.

Value

data.frame containing the segmentation.

Author(s)

Markus Riester

References

Olshen, A. B., Venkatraman, E. S., Lucito, R., Wigler, M. (2004). Circular binary segmentation for the analysis of array-based DNA copy number data. Biostatistics 5: 557-572.

Venkatraman, E. S., Olshen, A. B. (2007). A faster circular binary segmentation algorithm for the analysis of array CGH data. Bioinformatics 23: 657-63.

Olshen et al. (2011). Parent-specific copy number in paired tumor-normal studies using circular binary segmentation. Bioinformatics.

See Also

runAbsoluteCN

Examples


normal.coverage.file <- system.file("extdata", "example_normal_tiny.txt",
    package = "PureCN")
tumor.coverage.file <- system.file("extdata", "example_tumor_tiny.txt",
    package = "PureCN")
vcf.file <- system.file("extdata", "example.vcf.gz",
    package = "PureCN")

# The max.candidate.solutions, max.ploidy and test.purity parameters are set to
# non-default values to speed-up this example.  This is not a good idea for real
# samples.
 ret <-runAbsoluteCN(normal.coverage.file = normal.coverage.file,
     tumor.coverage.file = tumor.coverage.file, vcf.file = vcf.file,
     sampleid = "Sample1",  genome = "hg19",
     fun.segmentation = segmentationPSCBS, max.ploidy = 4,
     test.purity = seq(0.3, 0.7, by = 0.05), max.candidate.solutions = 1)


lima1/PureCN documentation built on Sept. 17, 2024, 5:48 a.m.