Title: | Comprehensive Pipeline for Analyzing and Visualizing Array-Based CGH Data |
---|---|
Description: | A comprehensive pipeline for analyzing and interactively visualizing genomic profiles generated through commercial or custom aCGH arrays. As inputs, rCGH supports Agilent dual-color Feature Extraction files (.txt), from 44 to 400K, Affymetrix SNP6.0 and cytoScanHD probeset.txt, cychp.txt, and cnchp.txt files exported from ChAS or Affymetrix Power Tools. rCGH also supports custom arrays, provided data complies with the expected format. This package takes over all the steps required for individual genomic profiles analysis, from reading files to profiles segmentation and gene annotations. This package also provides several visualization functions (static or interactive) which facilitate individual profiles interpretation. Input files can be in compressed format, e.g. .bz2 or .gz. |
Authors: | Frederic Commo [aut, cre] |
Maintainer: | Frederic Commo <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.37.0 |
Built: | 2024-12-06 06:15:49 UTC |
Source: | https://github.com/bioc/rCGH |
A comprehensive pipeline for analyzing and interactively visualizing genomic
profiles generated through commercial or custom aCGH arrays. As inputs, rCGH
supports Agilent dual-color Feature Extraction files (.txt), from 44 to
400K, Affymetrix SNP6.0 and cytoScanHD probeset.txt, cychp.txt, and
cnchp.txt files exported from ChAS or Affymetrix Power Tools.
rCGH also supports custom arrays, provided data is
in a suitable format. This package takes over all the steps required for
individual genomic profiles analysis, from reading files to segmenting
and annotating genes. This package provides several visualization functions
(static or interactive) which facilitate individual profiles interpretation.
Input files can be in compressed format, e.g. .bz2 or .gz.
Frederic Commo <[email protected]>
readAgilent
, readAffySNP6
,
readAffyCytoScan
, readGeneric
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Static visalizations plotDensity(cgh) multiplot(cgh) ## Not run: # Interactive visalizations view(cgh) ## End(Not run)
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Static visalizations plotDensity(cgh) multiplot(cgh) ## Not run: # Interactive visalizations view(cgh) ## End(Not run)
This function performs several preprocessing steps: local regressions
(loessFit) are used to correct cy3/cy5 and GC% bias, when Agilent dual-color
hybridization are used only.
In case of Affymetrix cychp.txt (or cnchp.txt) data are used, these steps have
been already done in ChAS or Affymetrix Power Tools.
## S4 method for signature 'rCGH' adjustSignal(object, Scale=TRUE, Cy=TRUE, GC=TRUE, Ref="cy3", suppOutliers=TRUE, nCores=NULL, verbose=TRUE)
## S4 method for signature 'rCGH' adjustSignal(object, Scale=TRUE, Cy=TRUE, GC=TRUE, Ref="cy3", suppOutliers=TRUE, nCores=NULL, verbose=TRUE)
object |
: An object of class |
Scale |
: logical. If |
Cy |
: logical. If |
GC |
: logical. If |
Ref |
: string. Which cyanine was used as the reference. Possible values are "cy3" (default) and "cy5". For Agilent dual-color hybridization only. |
suppOutliers |
: logical. If |
nCores |
: numeric. The number of cores to use in order to speed up the computation.
When |
verbose |
: logical. When |
When suppOutliers
is TRUE
(default), the Log2Ratios are splitted
with respect to chromosomes. The main regions within each chromosome are
identified using EM. Within each region r_i, x[r_i] values are redifined
according to the corresponding model parameters.
as:
Notice that this step may substantially increase the computation time.
An object of class "rCGH"
Frederic Commo
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) getParam(cgh)
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) getParam(cgh)
A dataset containing the Agilent probe Ids, and their GC content.
These data allow rCGH
to support Agilent aCGH arrays, from 44K to
400K. See the source
section.
This information is used to correct the GC% bias. For Agilent data only.
agilentDB
agilentDB
A data frame with 411056 rows and 2 columns:
ProbeID: Official Agilent probe ids.
GC: The GC content in each probe sequence, expressed as the GC fraction.
A dataset
Probe sequences are not used and have been removed after computing the GC fractions as GC = sum(G or C)/length(sequence), for each sequence in file.
Frederic Commo
These data derive from the official Agilent Sequence List,
SurePrint_G3_Human_CGH_Microarray 2x400K_021850_D_SequenceList_20111015.txt,
freely available at:
Agilent SureDesign
Access date: 3-2-2015
Notice that a User ID and Password are required to sign in.
"rCGH"
Accessor Functions Methods for extracting information from an object of
class "rCGH"
.
Each of the below methods are simply convenience functions which extract the
corresponding slots (as the name of each method suggests) from an object of
class "rCGH"
.
## S4 method for signature 'rCGH' getInfo(object, item = NULL) ## S4 method for signature 'rCGH' getCNset(object) ## S4 method for signature 'rCGH' getParam(object) ## S4 method for signature 'rCGH' getSegTable(object, minLen = NULL)
## S4 method for signature 'rCGH' getInfo(object, item = NULL) ## S4 method for signature 'rCGH' getCNset(object) ## S4 method for signature 'rCGH' getParam(object) ## S4 method for signature 'rCGH' getSegTable(object, minLen = NULL)
object |
: An object of class |
item |
: character. Can be one, or a vector of items.
When |
minLen |
: numeric. The mininal length for a segment, in Kb.
When |
getInfo(object, item = NULL): character.
getCNset(object): a data frame.
getParam(object): a list of parameters.
getSegTable(object, minLen = NULL): a data frame.
"rCGH"
getInfo(object, item = NULL): returns the values of the specified items, all the information otherwise.
getCNset(object): returns the full by-probe dataset.
getParam(object): returns the analysis parameters.
getSegTable(object, minLen = NULL): returns the segmentation table - one row per segment.
Frederic Commo
filePath <- system.file("extdata", "Agilent4x180K.txt.bz2", package = "rCGH") cgh <- readAgilent(filePath, sampleName = "Agilent4x180K", labName = "myLab") # Getting all the information getInfo(cgh) # Getting specific items getInfo(cgh, c("sampleName", "labName"))
filePath <- system.file("extdata", "Agilent4x180K.txt.bz2", package = "rCGH") cgh <- readAgilent(filePath, sampleName = "Agilent4x180K", labName = "myLab") # Getting all the information getInfo(cgh) # Getting specific items getInfo(cgh, c("sampleName", "labName"))
This function creates a by-gene table by listing all the genes contained in each
of the segments in the segmentation table.
Gene annotations (symbol, location,...), segmented Log2Ratios, and segment
legnth are reported in the final table.
A supplementary score is the relativeLog
: the magnitude, in Log2, from
the closest centromere.
byGeneTable(segTable, symbol = NULL, genome = c("hg19", "hg18", "hg38"), columns = NA, verbose = TRUE)
byGeneTable(segTable, symbol = NULL, genome = c("hg19", "hg18", "hg38"), columns = NA, verbose = TRUE)
segTable |
: data frame. A segmentation table exported from an object of class
|
symbol |
: character. A valid HUGO symbol. When |
genome |
: string. The genome build to use. Supported genomes are |
columns |
: string. what supplementary genes annotations to export.
Allowed annotations are those supported by the
|
verbose |
: logical. When |
For gene annotations, Hg19/GRCh37 annotations downloaded from NCBI are considered.
An object of class "rCGH"
Frederic Commo
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) st <- getSegTable(cgh) bygene <- byGeneTable(st) head(bygene)
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) st <- getSegTable(cgh) bygene <- byGeneTable(st) head(bygene)
This function analyses the Log2Ratios as a mixture of several gaussian
populations, using an Expectation-Maximization algorithm (EM).
The peakThresh
argument specifies what proportion of the main density
peak is allowed for choosing a neutral 2-copies population. The mean of the
chosen population is used for centralizing the profile.
See Mclust
.
## S4 method for signature 'rCGH' EMnormalize(object, G = 2:6, priorScale = 5, peakThresh = 0.5, mergeVal = 0.1, Title = NA, verbose = TRUE)
## S4 method for signature 'rCGH' EMnormalize(object, G = 2:6, priorScale = 5, peakThresh = 0.5, mergeVal = 0.1, Title = NA, verbose = TRUE)
object |
: An object of class |
G |
: numeric. The number of groups to test during the gaussian mixture estimation. Default is from 2 to 6. |
priorScale |
: numeric. A scale value passed to |
peakThresh |
: numeric. The proportion of the highest peak to consider as a peak selection threshold. Default is 0.5. |
mergeVal |
: numeric. Populations with means closer than |
Title |
: character string. A title for the density plot. If |
verbose |
: logical. When |
Depending on peakThresh
, the mean of the highest density, or a lower
value, is chosen for centering the genomic profile. To do so, L2R are modeled
for each segment s_i, with respect to n_i (the number of probes included in
segment i), mu_i and sd_i. The mixture of L2Rs is then analysed as a mixture of
gaussian populations.
When a peakThresh
value is specified, heights of density peaks are
compared: the lowest peak mean among the peaks respecting the criteria:
peakHeight > max(peaks)*peakThresh
, is chosen for centralizing the data.
See References
An object of same class as the input.
Frederic Commo
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) getParam(cgh)
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) getParam(cgh)
A data set containing lengths and centromere locations for each of the 24 chromosomes, according to Hg18.
hg18
hg18
A data set with 24 rows and 5 columns:
chrom: chromosome number.
length: chromosome length.
centromerStart: centromere start position.
centromerEnd: centromere start position.
cumlen: cumulative length (where the previous chromosome ends).
a data set.
Frederic Commo
These data derived from the Hg18 gap UCSC table, freely available at:
UCSC
Access date: 10-10-2015
Within the browser, select:
group: All Tables
database: hg18
table: gap
# For users convenience, we provide a prebuilt dataset # containing the Hg18 chr lengths, and centromeres location. hg18 # The same dataset can be obtained as follow: ## Not run: library(BSgenome) library(rtracklayer) getChrLength <- function(genome){ genome <- sprintf("BSgenome.Hsapiens.UCSC. g <- getBSgenome(genome, masked=FALSE) data.frame(chrom=1:24, length=seqlengths(g)[1:24]) } .chrAsNum <- function(tbl){ tbl$chrom <- gsub("chr", "", tbl$chrom) tbl$chrom[tbl$chrom=="X"] <- 23 tbl$chrom[tbl$chrom=="Y"] <- 24 tbl$chrom <- as.numeric(tbl$chrom) tbl[order(tbl$chrom),] } getCentromeres <- function(genome){ mySession <- try(browserSession("UCSC"), silent=TRUE) # In case it fails, use another mirror if(inherits(mySession, "try-error")) mySession <- browserSession("UCSC", url="http://genome-euro.ucsc.edu/cgi-bin/") genome(mySession) <- genome obj <- ucscTableQuery(mySession, table="gap") tbl <- getTable(obj) tbl <- tbl[tbl$type=="centromere", c("chrom", "chromStart", "chromEnd")] colnames(tbl)[2:3] <- c("centromerStart", "centromerEnd") .chrAsNum(tbl) } makeHg <- function(genome){ chrL <- getChrLength(genome) ctm <- getCentromeres(genome) tbl <- merge(chrL, ctm, by="chrom") cumlen <- c(0, cumsum(as.numeric(tbl$length))[-nrow(tbl)]) cbind.data.frame(tbl, cumlen=cumlen) } hg18 <- makeHg("hg18") hg18 ## End(Not run)
# For users convenience, we provide a prebuilt dataset # containing the Hg18 chr lengths, and centromeres location. hg18 # The same dataset can be obtained as follow: ## Not run: library(BSgenome) library(rtracklayer) getChrLength <- function(genome){ genome <- sprintf("BSgenome.Hsapiens.UCSC. g <- getBSgenome(genome, masked=FALSE) data.frame(chrom=1:24, length=seqlengths(g)[1:24]) } .chrAsNum <- function(tbl){ tbl$chrom <- gsub("chr", "", tbl$chrom) tbl$chrom[tbl$chrom=="X"] <- 23 tbl$chrom[tbl$chrom=="Y"] <- 24 tbl$chrom <- as.numeric(tbl$chrom) tbl[order(tbl$chrom),] } getCentromeres <- function(genome){ mySession <- try(browserSession("UCSC"), silent=TRUE) # In case it fails, use another mirror if(inherits(mySession, "try-error")) mySession <- browserSession("UCSC", url="http://genome-euro.ucsc.edu/cgi-bin/") genome(mySession) <- genome obj <- ucscTableQuery(mySession, table="gap") tbl <- getTable(obj) tbl <- tbl[tbl$type=="centromere", c("chrom", "chromStart", "chromEnd")] colnames(tbl)[2:3] <- c("centromerStart", "centromerEnd") .chrAsNum(tbl) } makeHg <- function(genome){ chrL <- getChrLength(genome) ctm <- getCentromeres(genome) tbl <- merge(chrL, ctm, by="chrom") cumlen <- c(0, cumsum(as.numeric(tbl$length))[-nrow(tbl)]) cbind.data.frame(tbl, cumlen=cumlen) } hg18 <- makeHg("hg18") hg18 ## End(Not run)
A data set containing lengths and centromere locations for each of the 24 chromosomes, according to Hg19.
hg19
hg19
A data set with 24 rows and 5 columns:
chrom: chromosome number.
length: chromosome length.
centromerStart: centromere start position.
centromerEnd: centromere start position.
cumlen: cumulative length (where the previous chromosome ends).
a data set.
Frederic Commo
These data derived from the Hg19 gap UCSC table, freely available at:
UCSC
Access date: 1-31-2014
Within the browser, select:
group: All Tables
database: hg19
table: gap
# For users convenience, we provide a prebuilt dataset # containing the Hg19 chr lengths, and centromeres location. hg19 # The same dataset can be obtained as follow: ## Not run: library(BSgenome) library(rtracklayer) getChrLength <- function(genome = "BSgenome.Hsapiens.UCSC.hg19"){ g <- getBSgenome(genome, masked=FALSE) data.frame(chrom=1:24, length=seqlengths(g)[1:24]) } .chrAsNum <- function(tbl){ tbl$chrom <- gsub("chr", "", tbl$chrom) tbl$chrom[tbl$chrom=="X"] <- 23 tbl$chrom[tbl$chrom=="Y"] <- 24 tbl$chrom <- as.numeric(tbl$chrom) tbl[order(tbl$chrom),] } getCentromeres <- function( genome="hg19" ){ mySession <- try(browserSession("UCSC"), silent=TRUE) # In case of failure, try another mirror if(inherits(mySession, "try-error")) mySession <- browserSession("UCSC", url="http://genome-euro.ucsc.edu/cgi-bin/") genome(mySession) <- genome obj <- ucscTableQuery(mySession, table="gap") tbl <- getTable(obj) tbl <- tbl[tbl$type=="centromere", c("chrom", "chromStart", "chromEnd")] colnames(tbl)[2:3] <- c("centromerStart", "centromerEnd") .chrAsNum(tbl) } makeHg19 <- function(){ tbl <- merge(getChrLength(), getCentromeres(), by="chrom") cumlen <- c(0, cumsum(as.numeric(tbl$length))[-nrow(tbl)]) cbind.data.frame(tbl, cumlen=cumlen) } hg19 <- makeHg19() hg19 ## End(Not run)
# For users convenience, we provide a prebuilt dataset # containing the Hg19 chr lengths, and centromeres location. hg19 # The same dataset can be obtained as follow: ## Not run: library(BSgenome) library(rtracklayer) getChrLength <- function(genome = "BSgenome.Hsapiens.UCSC.hg19"){ g <- getBSgenome(genome, masked=FALSE) data.frame(chrom=1:24, length=seqlengths(g)[1:24]) } .chrAsNum <- function(tbl){ tbl$chrom <- gsub("chr", "", tbl$chrom) tbl$chrom[tbl$chrom=="X"] <- 23 tbl$chrom[tbl$chrom=="Y"] <- 24 tbl$chrom <- as.numeric(tbl$chrom) tbl[order(tbl$chrom),] } getCentromeres <- function( genome="hg19" ){ mySession <- try(browserSession("UCSC"), silent=TRUE) # In case of failure, try another mirror if(inherits(mySession, "try-error")) mySession <- browserSession("UCSC", url="http://genome-euro.ucsc.edu/cgi-bin/") genome(mySession) <- genome obj <- ucscTableQuery(mySession, table="gap") tbl <- getTable(obj) tbl <- tbl[tbl$type=="centromere", c("chrom", "chromStart", "chromEnd")] colnames(tbl)[2:3] <- c("centromerStart", "centromerEnd") .chrAsNum(tbl) } makeHg19 <- function(){ tbl <- merge(getChrLength(), getCentromeres(), by="chrom") cumlen <- c(0, cumsum(as.numeric(tbl$length))[-nrow(tbl)]) cbind.data.frame(tbl, cumlen=cumlen) } hg19 <- makeHg19() hg19 ## End(Not run)
A data set containing lengths and centromere locations for each of the 24 chromosomes, according to Hg38.
hg38
hg38
A data set with 24 rows and 5 columns:
chrom: chromosome number.
length: chromosome length.
centromerStart: centromere start position.
centromerEnd: centromere start position.
cumlen: cumulative length (where the previous chromosome ends).
a data set.
Frederic Commo
These data derived from the Hg38 gap UCSC table, freely available at:
UCSC
Access date: 10-10-2015
Within the browser, select:
group: All Tables
database: hg38
table: gap
# For users convenience, we provide a prebuilt dataset # containing the Hg38 chr lengths, and centromeres location. hg38 # The same dataset can be obtained as follow: ## Not run: library(BSgenome) library(rtracklayer) getChrLength <- function(genome){ genome <- sprintf("BSgenome.Hsapiens.UCSC. g <- getBSgenome(genome, masked=FALSE) data.frame(chrom=1:24, length=seqlengths(g)[1:24]) } .chrAsNum <- function(tbl){ tbl$chrom <- gsub("chr", "", tbl$chrom) tbl$chrom[tbl$chrom=="X"] <- 23 tbl$chrom[tbl$chrom=="Y"] <- 24 tbl$chrom <- as.numeric(tbl$chrom) tbl[order(tbl$chrom),] } getCentromeres <- function(genome){ mySession <- try(browserSession("UCSC"), silent=TRUE) # In case it fails, use another mirror if(inherits(mySession, "try-error")) mySession <- browserSession("UCSC", url="http://genome-euro.ucsc.edu/cgi-bin/") genome(mySession) <- genome obj <- ucscTableQuery(mySession, table="gap") tbl <- getTable(obj) if(!"centromere" return(NULL) tbl <- tbl[tbl$type=="centromere", c("chrom", "chromStart", "chromEnd")] colnames(tbl)[2:3] <- c("centromerStart", "centromerEnd") .chrAsNum(tbl) } makeHg <- function(genome){ chrL <- getChrLength(genome) ctm <- getCentromeres(genome) # Notice that, in case of Hg38, centromeres locations are in Hg19. if(is.null(ctm)) ctm <- getCentromeres("hg19") tbl <- merge(chrL, ctm, by="chrom") cumlen <- c(0, cumsum(as.numeric(tbl$length))[-nrow(tbl)]) cbind.data.frame(tbl, cumlen=cumlen) } hg38 <- makeHg("hg38") hg38 ## End(Not run)
# For users convenience, we provide a prebuilt dataset # containing the Hg38 chr lengths, and centromeres location. hg38 # The same dataset can be obtained as follow: ## Not run: library(BSgenome) library(rtracklayer) getChrLength <- function(genome){ genome <- sprintf("BSgenome.Hsapiens.UCSC. g <- getBSgenome(genome, masked=FALSE) data.frame(chrom=1:24, length=seqlengths(g)[1:24]) } .chrAsNum <- function(tbl){ tbl$chrom <- gsub("chr", "", tbl$chrom) tbl$chrom[tbl$chrom=="X"] <- 23 tbl$chrom[tbl$chrom=="Y"] <- 24 tbl$chrom <- as.numeric(tbl$chrom) tbl[order(tbl$chrom),] } getCentromeres <- function(genome){ mySession <- try(browserSession("UCSC"), silent=TRUE) # In case it fails, use another mirror if(inherits(mySession, "try-error")) mySession <- browserSession("UCSC", url="http://genome-euro.ucsc.edu/cgi-bin/") genome(mySession) <- genome obj <- ucscTableQuery(mySession, table="gap") tbl <- getTable(obj) if(!"centromere" return(NULL) tbl <- tbl[tbl$type=="centromere", c("chrom", "chromStart", "chromEnd")] colnames(tbl)[2:3] <- c("centromerStart", "centromerEnd") .chrAsNum(tbl) } makeHg <- function(genome){ chrL <- getChrLength(genome) ctm <- getCentromeres(genome) # Notice that, in case of Hg38, centromeres locations are in Hg19. if(is.null(ctm)) ctm <- getCentromeres("hg19") tbl <- merge(chrL, ctm, by="chrom") cumlen <- c(0, cumsum(as.numeric(tbl$length))[-nrow(tbl)]) cbind.data.frame(tbl, cumlen=cumlen) } hg38 <- makeHg("hg38") hg38 ## End(Not run)
This function display a static view of the genomic profile and the allelic
difference stored in an object of class "rCGH"
.
If no allelic difference is available, the genomic profile only is displayed.
## S4 method for signature 'rCGH' multiplot(object, symbol=NULL, gain=.5, loss=(-.5), minLen=10, pCol = "grey50", GLcol = c("blue", "red3"), L=matrix(seq(1, 12)), p=c(1/2, 1/4, 1/4), Title=NULL, ylim=NULL)
## S4 method for signature 'rCGH' multiplot(object, symbol=NULL, gain=.5, loss=(-.5), minLen=10, pCol = "grey50", GLcol = c("blue", "red3"), L=matrix(seq(1, 12)), p=c(1/2, 1/4, 1/4), Title=NULL, ylim=NULL)
object |
: An object of class |
symbol |
: character. A valid HUGO symbol (case insensitive). |
gain |
: numeric. A gain threshold value (in |
loss |
: numeric. A loss threshold value (in |
minLen |
: numeric. The mininal length for a segment, expressed in Kb.
When |
pCol |
: string. The probe points color. Default is |
GLcol |
: vector. A vector of 2 colors: the gained and lost segments colors,
respectively. Default is |
L |
: matrix. A matrix defining the layout. Default is 12 lines. |
p |
: numeric. The proportion of each plot within the plot window.
Default is |
Title |
: character string. A title for the plot. If |
ylim |
: numeric. A vector of two values specifying the y-axis range.
See |
None.
If no allelic difference is available, the genomic profile only is displayed.
Frederic Commo
plotDensity
, plotProfile
,
plotLOH
, view
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Static visalizations multiplot(cgh, symbol = "erbb2")
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Static visalizations multiplot(cgh, symbol = "erbb2")
This function display the distribution of the Log2Ratios, as well as how the
"EMnormalize"
step estimates the mixture of gaussian
populations, and choose a centralization value.
## S4 method for signature 'rCGH' plotDensity(object, breaks=NULL, Title=NULL,...)
## S4 method for signature 'rCGH' plotDensity(object, breaks=NULL, Title=NULL,...)
object |
: An object of class |
breaks |
: The number of breaks to use. See |
Title |
: character string. A title for the density plot. If |
... |
: Other graphical parameters supported by |
None.
Frederic Commo
plotProfile
, plotLOH
, multiplot
,
view
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) plotDensity(cgh)
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) plotDensity(cgh)
This function display a static view of the allele differences, when available.
## S4 method for signature 'rCGH' plotLOH(object, Title=NULL)
## S4 method for signature 'rCGH' plotLOH(object, Title=NULL)
object |
: An object of class |
Title |
: character string. A title for the density plot. If |
None.
Frederic Commo
plotDensity
, plotProfile
, multiplot
,
view
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Static visalizations plotLOH(cgh)
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Static visalizations plotLOH(cgh)
This function display a static view of the genomic profile stored in an
object of class "rCGH"
.
## S4 method for signature 'rCGH' plotProfile(object, showCopy=FALSE, symbol=NULL, gain=.5, loss=(-.5), minLen = 10, pCol = "grey50", GLcol = c("blue", "red3"), Title=NULL, ylim=NULL)
## S4 method for signature 'rCGH' plotProfile(object, showCopy=FALSE, symbol=NULL, gain=.5, loss=(-.5), minLen = 10, pCol = "grey50", GLcol = c("blue", "red3"), Title=NULL, ylim=NULL)
object |
: An object of class |
showCopy |
: logical. To show the estimated copy numbers instead of the Log2Ratios.
default is |
symbol |
: character. A valid HUGO symbol (case insensitive). |
gain |
: numeric. A gain threshold value (in |
loss |
: numeric. A loss threshold value (in |
minLen |
: numeric. The mininal length for a segment, in Kb.
When |
pCol |
: string. The probe points color. DEfault is |
GLcol |
: vector. A vector of 2 colors: the gained and lost segments colors,
respectively. Default is |
Title |
: character string. A title for the density plot. If |
ylim |
: numeric. A vector of two values specifying a range for the y-axis.
If |
None.
Frederic Commo
plotDensity
, plotLOH
, multiplot
,
view
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Static visalization using Log2Ratios plotProfile(cgh, symbol = "erbb2") # Static visalization using estimated copy numbers plotProfile(cgh, showCopy = TRUE, symbol = "erbb2")
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Static visalization using Log2Ratios plotProfile(cgh, symbol = "erbb2") # Static visalization using estimated copy numbers plotProfile(cgh, showCopy = TRUE, symbol = "erbb2")
"rCGH-Agilent"
An instance of class "rCGH-Agilent"
, which inherits from
the superclass "rCGH"
.
Slots described below are used to store sample information, analysis
parameters, and segmentation results. All are accessible through specific
"Accessors"
functions.
Objects can be created by calls of the form new("rCGH-Agilent", ...)
.
info
:Object of class "character"
: where sample information can be
stored. See getInfo
and setInfo
.
cnSet
:Object of class "data.frame"
: the full data set.
See getCNset
.
param
:Object of class "list"
: the analysis parameters stored for
traceability. getParam
.
segTable
:Object of class "data.frame"
: the segmentation table.
getSegTable
.
Class "rCGH"
, directly.
No methods defined with class "rCGH-Agilent" in the signature.
Frederic Commo
"rCGH"
, "rCGH-cytoScan"
,
"rCGH-SNP6"
, "rCGH-generic"
showClass("rCGH-Agilent")
showClass("rCGH-Agilent")
"rCGH"
Class "rCGH"
is a superclass living on top of
"rCGH-Agilent"
, "rCGH-SNP6"
,
"rCGH-cytoScan"
, "rCGH-oncoScan"
,
and "rCGH-generic"
.
These objects inherit most of the properties of the superclass, and allow
specific parameterizations used during the analysis process.
Objects are created by platform-specific read
functions:
"readAgilent"
, "readAffySNP6"
,
and "readAffyCytoScan"
,
each corresponding to their matched file format.
A supplementary "readGeneric"
allows the user to create a
"rCGH"
object from custom arrays.
Slots described below are used to store sample information and analysis
parameters, as well as segmentation results. All are accessible through specific
"Accessors"
functions.
Objects can be created by calls of the form new("rCGH", ...)
.
Slots content are updated at each different analysis step, and are accessible
through specific get
functions.
info
:Object of class "character"
: where sample information can be
stored. See "getInfo"
and "setInfo"
.
cnSet
:Object of class "data.frame"
: the full data set.
See "getCNset"
.
param
:Object of class "list"
: the analysis parameters stored for
traceability. See "getParam"
.
segTable
:Object of class "data.frame"
: the segmentation table.
See "getSegTable"
.
signature(object = "rCGH")
: ...
Frederic Commo
"rCGH-Agilent"
, "rCGH-SNP6"
,
"rCGH-cytoScan"
, "rCGH-generic"
showClass("rCGH")
showClass("rCGH")
"rCGH-cytoScan"
An instance of class "rCGH-cytoScan"
, which inherits from
the superclass "rCGH"
.
Slots described below are used to store sample information, analysis
parameters, and segmentation results. All are accessible through specific
"Accessors"
functions.
Objects can be created by calls of the form new("rCGH-cytoScan", ...)
.
info
:Object of class "character"
: where sample information can be
stored. See getInfo
and setInfo
.
cnSet
:Object of class "data.frame"
: the full data set.
See getCNset
.
param
:Object of class "list"
: the analysis parameters stored for
traceability. getParam
.
segTable
:Object of class "data.frame"
: the segmentation table.
getSegTable
.
Class "rCGH"
, directly.
No methods defined with class "rCGH-cytoScan" in the signature.
Frederic Commo
"rCGH"
, "rCGH-Agilent"
,
"rCGH-SNP6"
, "rCGH-generic"
showClass("rCGH-cytoScan")
showClass("rCGH-cytoScan")
"rCGH-generic"
An instance of class "rCGH-generic"
, which inherits from
the superclass "rCGH"
.
Slots described below are used to store sample information, analysis
parameters, and segmentation results. All are accessible through specific
"Accessors"
functions.
Objects can be created by calls of the form new("rCGH-generic", ...)
.
info
:Object of class "character"
: where sample information can be
stored. See getInfo
and setInfo
.
cnSet
:Object of class "data.frame"
: the full data set.
See getCNset
.
param
:Object of class "list"
: the analysis parameters stored for
traceability. getParam
.
segTable
:Object of class "data.frame"
: the segmentation table.
getSegTable
.
Class "rCGH"
, directly.
No methods defined with class "rCGH-generic" in the signature.
Frederic Commo
"rCGH"
, "rCGH-Agilent"
,
"rCGH-SNP6"
, "rCGH-cytoScan"
showClass("rCGH-generic")
showClass("rCGH-generic")
"rCGH-oncoScan"
An instance of class "rCGH-oncoScan"
, which inherits from
the superclass "rCGH"
.
Slots described below are used to store sample information, analysis
parameters, and segmentation results. All are accessible through specific
"Accessors"
functions.
Objects can be created by calls of the form new("rCGH-oncoScan", ...)
.
info
:Object of class "character"
: where sample information can be
stored. See getInfo
and setInfo
.
cnSet
:Object of class "data.frame"
: the full data set.
See getCNset
.
param
:Object of class "list"
: the analysis parameters stored for
traceability. getParam
.
segTable
:Object of class "data.frame"
: the segmentation table.
getSegTable
.
Class "rCGH"
, directly.
No methods defined with class "rCGH-oncoScan" in the signature.
Frederic Commo
"rCGH"
, "rCGH-Agilent"
,
"rCGH-SNP6"
, "rCGH-cytoScan"
,
"rCGH-generic"
showClass("rCGH-oncoScan")
showClass("rCGH-oncoScan")
"rCGH-SNP6"
An instance of class "rCGH-SNP6"
, which inherits from
the superclass "rCGH"
.
Slots described below are used to store sample information, analysis
parameters, and segmentation results. All are accessible through specific
"Accessors"
functions.
Objects can be created by calls of the form new("rCGH-SNP6", ...)
.
info
:Object of class "character"
: where sample information can be
stored. See getInfo
and setInfo
.
cnSet
:Object of class "data.frame"
: the full data set.
See getCNset
.
param
:Object of class "list"
: the analysis parameters stored for
traceability. getParam
.
segTable
:Object of class "data.frame"
: the segmentation table.
getSegTable
.
Class "rCGH"
, directly.
No methods defined with class "rCGH-SNP6" in the signature.
Frederic Commo
"rCGH"
, "rCGH-Agilent"
,
"rCGH-cytoScan"
showClass("rCGH-SNP6")
showClass("rCGH-SNP6")
"rCGH-cytoScan"
ConstructorA constructor function which takes an Affymetrix cytoScanHD cychp.txt
(or cnchp.txt) file as input, possibly in a compressed format (.bz2 or .gz).
These files are exported from Chromosome Analysis Suite (ChAS) or Affymetrix
Power Tools (see References
section).
readAffyCytoScan(filePath, sampleName=NA, labName=NA, useProbes=c("snp", "cn", "all"), genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose=TRUE)
readAffyCytoScan(filePath, sampleName=NA, labName=NA, useProbes=c("snp", "cn", "all"), genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose=TRUE)
filePath |
: string. A path to an Affymetrix cytoScanHD cychp.txt (or cnchp.txt) file. |
sampleName |
: string. A sample Id. Optional. |
labName |
: string. A lab Id. Optional. |
useProbes |
: character. What probes to consider. Possible choices are SNP probes only ("snp", default), CN probes only ("cn"), or all the probes ("all"). |
genome |
: string. The genome build to use. Supported genomes are |
ploidy |
: numeric. A priori plody value, when known, to adjust the estimation of copy numbers. Default is 2. |
verbose |
: logical. When |
When available in the file preambule, several array information will be stored
in Object@info
: scanning date, grid version,...
Any other useful item can be stored using setInfo
.
An object of class "rCGH"
Frederic Commo
readAgilent
, readAffySNP6
,
readGeneric
, readAffyOncoScan
,
setInfo
, getInfo
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh
"rCGH-oncoScan"
ConstructorA constructor function which takes an Affymetrix oncoScan tabulated
file as input, possibly in a compressed format (.bz2 or .gz).
This can be either a 'ProbeSets,CopyNumber.tsv'
alone, or merged with its
corresponding 'ProbeSets,AllelicData.tsv'
file. See the details
section.
readAffyOncoScan(filePath, sampleName=NA, labName=NA, genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose=TRUE)
readAffyOncoScan(filePath, sampleName=NA, labName=NA, genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose=TRUE)
filePath |
: string. A path to an Affymetrix .tsv file. See |
sampleName |
: string. A sample Id. Optional. |
labName |
: string. A lab Id. Optional. |
genome |
: string. The genome build to use. Supported genomes are |
ploidy |
: numeric. A priori plody value, when known, to adjust the estimation of copy numbers. Default is 2. |
verbose |
: logical. When |
The Affymetrix Power Tools apt-copynumber-onco-ssa
script produces 2
files: 'ProbeSets,CopyNumber.tsv'
and 'ProbeSets,AllelicData.tsv'
.
Merging these 2 files may produce a unique file containing both probes
Log2Ratio
and AllelicDifference
.
An object of class "rCGH"
Frederic Commo
readAgilent
, readAffySNP6
,
readGeneric
, readAffyCytoScan
,
setInfo
, getInfo
# Just a toy file filePath <- system.file("extdata", "oncoscan.tsv.bz2", package = "rCGH") cgh <- readAffyOncoScan(filePath, sampleName = "AffyOncoScan") cgh
# Just a toy file filePath <- system.file("extdata", "oncoscan.tsv.bz2", package = "rCGH") cgh <- readAffyOncoScan(filePath, sampleName = "AffyOncoScan") cgh
"rCGH-SNP6"
ConstructorA constructor function which takes an Affymetrix SNP6 cychp.txt (or cnchp.txt)
file as input, possibly in a compressed format (.bz2 or .gz).
These files are exported from Chromosome Analysis Suite (ChAS) or Affymetrix
Power Tools (APT) (see the References
section).
readAffySNP6(filePath, sampleName = NA, labName = NA, useProbes=c("snp", "cn", "all"), genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose = TRUE)
readAffySNP6(filePath, sampleName = NA, labName = NA, useProbes=c("snp", "cn", "all"), genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose = TRUE)
filePath |
: string. A path to an Affymetrix SNP6 cychp.txt (or cnchp.txt) file. |
sampleName |
: string. A sample Id. Optional. |
labName |
: string. A lab Id. Optional. |
useProbes |
: character. What probes to consider. Possible choices are SNP probes only ("snp", default), CN probes only ("cn"), or all the probes ("all"). |
genome |
: string. The genome build to use. Supported genomes are |
ploidy |
: numeric. A priori ploidy value, when known, to adjust the estimation of copy numbers. Default is 2. |
verbose |
: logical. When |
When available in the file preambule, several array information will be stored
in Object@info
: scanning date, grid version,...
Any other useful item can be stored using setInfo
.
An object of class "rCGH"
Frederic Commo
readAgilent
, readAffyCytoScan
,
readGeneric
, readAffyOncoScan
,
setInfo
, getInfo
filePath <- system.file("extdata", "Affy_snp6_cnchp.txt.bz2", package = "rCGH") cgh <- readAffySNP6(filePath, sampleName = "AffySNP6") cgh
filePath <- system.file("extdata", "Affy_snp6_cnchp.txt.bz2", package = "rCGH") cgh <- readAffySNP6(filePath, sampleName = "AffySNP6") cgh
"rCGH-Agilent"
Constructor.
A constructor function taking as input an Agilent FE .txt file, exported from
Feature Extraction, possibly in a compressed format (.bz2 or .gz).
Agilent from 44 to 400K are supported.
readAgilent(filePath, sampleName = NA, labName = NA, supFlags = TRUE, genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose = TRUE)
readAgilent(filePath, sampleName = NA, labName = NA, supFlags = TRUE, genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose = TRUE)
filePath |
: string. A path to an Agilent FE (.txt) file. |
sampleName |
: string. A sample Id. Optional. |
labName |
: string. A lab Id. Optional. |
supFlags |
: should the flagged probes be suppressed. Default is |
genome |
: string. The genome build to use. Supported genomes are |
ploidy |
: numeric. A priori ploidy value, when known, to adjust the estimation of copy numbers. Default is 2. |
verbose |
: logical. if |
When available in the file preambule, several array information will be stored
in Object@info
: scanning date, grid version,...
Any other useful item can be stored using setInfo
.
An object of class "rCGH"
Frederic Commo
readAffyCytoScan
, readAffySNP6
,
readGeneric
, readAffyOncoScan
,
setInfo
, getInfo
filePath <- system.file("extdata", "Agilent4x180K.txt.bz2", package = "rCGH") cgh <- readAgilent(filePath, sampleName = "Agilent4x180K", labName = "myLab") cgh
filePath <- system.file("extdata", "Agilent4x180K.txt.bz2", package = "rCGH") cgh <- readAgilent(filePath, sampleName = "Agilent4x180K", labName = "myLab") cgh
"rCGH-generic"
ConstructorA constructor function which takes a tabulated .txt file as input, possibly in a
compressed format (.bz2 or .gz).
Notice that precise column names are mandatory, see the details
section.
readGeneric(filePath, sampleName=NA, labName=NA, genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose=TRUE)
readGeneric(filePath, sampleName=NA, labName=NA, genome = c("hg19", "hg18", "hg38"), ploidy = 2, verbose=TRUE)
filePath |
: string. A path to an Generic .txt file. |
sampleName |
: string. A sample Id. Optional. |
labName |
: string. A lab Id. Optional. |
genome |
: string. The genome build to use. Supported genomes are |
ploidy |
: numeric. A priori ploidy value, when known, to adjust the estimation of copy numbers. Default is 2. |
verbose |
: logical. When |
This generic constructor does not expect any preamble. Mandatory columns are:
ProbeName
:Character strings. Typicaly the probe ids.
ChrNum
:numeric. The chromosome numbers. In case Chr X and Y are used and named
as "X"
and "Y"
, these notations will be converted into
23 and 24, respectively.
ChrStart
:numeric. The chromosomal probes locations.
Log2Ratio
:numeric. The corresponding Log2Ratios.
An object of class "rCGH"
Frederic Commo
readAgilent
, readAffySNP6
,
readAffyCytoScan
, readAffyOncoScan
,
setInfo
, getInfo
filePath <- system.file("extdata", "generic.txt.bz2", package = "rCGH") cgh <- readGeneric(filePath, sampleName = "demo") cgh
filePath <- system.file("extdata", "generic.txt.bz2", package = "rCGH") cgh <- readGeneric(filePath, sampleName = "demo") cgh
This function allows the user to recenter a genomic profile stored in
an object of class "rCGH"
.
Peaks are indexed from 1 to k, from left to right, as they appear on the
plotDensity
after the EMnormalize
step.
## S4 replacement method for signature 'rCGH' recenter(object) <- value
## S4 replacement method for signature 'rCGH' recenter(object) <- value
object |
: An object of class |
value |
: numeric. What peak number to choose to recenter the genomic profile. |
An object of class "rCGH"
When a profile is recentered, the stored workflow parameters are updated.
see getParam
.
Frederic Commo
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Default peak choice center the profile on the 1st peak plotDensity(cgh) # Recentering on the 2nd density peak recenter(cgh) <- 2 plotDensity(cgh)
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) # Default peak choice center the profile on the 1st peak plotDensity(cgh) # Recentering on the 2nd density peak recenter(cgh) <- 2 plotDensity(cgh)
A function for performing the Log2Ratio segmentation on an object of class
"rCGH"
. See the details
section below.
## S4 method for signature 'rCGH' segmentCGH(object, Smooth=TRUE, UndoSD = NULL, minLen = 10, nCores=NULL, verbose = TRUE)
## S4 method for signature 'rCGH' segmentCGH(object, Smooth=TRUE, UndoSD = NULL, minLen = 10, nCores=NULL, verbose = TRUE)
object |
: An object of class |
Smooth |
: logical. Should the LRR be smoothed before being segmented.
See |
UndoSD |
: numeric. When not specified (default is |
minLen |
: numeric. The minimal length for a segment, in Kb. Shorter segments will be merged to the closest adjacent one. Default value is 10(Kb). |
nCores |
: numeric. The number of cores to use in order to speed up the computation.
When |
verbose |
: logical. if |
This function is a wrapper for the DNAcopy
, CNA
and segment
functions, which allows parallelization and
data-driven parameterization.
In addition to the usual DNAcopy
output, the segmentation table contains
the probes Log2Ratio standard deviation for each segment, as well as there
length, in Kb.
An object of class "rCGH"
Frederic Commo
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) st <- getSegTable(cgh) head(st)
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) st <- getSegTable(cgh) head(st)
"rCGH"
This function allows the user to store any type of supplementaty information in
an object of class "rCGH"
.
## S4 replacement method for signature 'rCGH' setInfo(object, item = NULL) <- value
## S4 replacement method for signature 'rCGH' setInfo(object, item = NULL) <- value
object |
: An object of class |
item |
: character. An item name to store. Default is |
value |
: any. A value to store. |
An object of class "rCGH"
When either item
or value
is NULL
, an error is returned.
Frederic Commo
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") # When supplementary information is added, # numerical, logical, or strings are supported setInfo(cgh, "someItem1") <- 35 setInfo(cgh, "someItem2") <- TRUE setInfo(cgh, "someItem3") <- "someComment" getInfo(cgh) # or to get back specific items getInfo(cgh, c("someItem1", "someItem3"))
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") # When supplementary information is added, # numerical, logical, or strings are supported setInfo(cgh, "someItem1") <- 35 setInfo(cgh, "someItem2") <- TRUE setInfo(cgh, "someItem3") <- "someComment" getInfo(cgh) # or to get back specific items getInfo(cgh, c("someItem1", "someItem3"))
show
"rCGH"
A method for visualizing an object of class "rCGH"
signature(object = "rCGH")
Frederic Commo
This function is build on top of shiny
, and provides an
iteractive way for visualizing a genomic profile, and exploring the list of
genes.
From a command panel, the user can interact with the graph in different
ways. See details
.
## S4 method for signature 'rCGH' view(object, browser = TRUE, ...)
## S4 method for signature 'rCGH' view(object, browser = TRUE, ...)
object |
: An object of class |
browser |
: logical. When |
... |
: Optional parameters used by |
The left command panel allows the user several actions:
displaying a specific gene by calling its HUGO symbol.
showing all or one unique chromosome.
merging segments shorter than a specified value, in Kb.
recentering the entire profile.
rescaling the y-axis.
specifying the Log2Ratio cut offs for defining gains and losses.
specifying a segment lenght cut off, in Mb.
exporting the genomic plot.
exporting the genes list.
Some actions, such as showing one unique chromosome or specifying cut offs (gain, loss, segment length), automatically update the gene table available in the "Genes table" tab.
None.
Frederic Commo
plotProfile
, plotLOH
, multiplot
,
runApp
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) ## Not run: # Interactive visalizations view(cgh) ## End(Not run)
filePath <- system.file("extdata", "Affy_cytoScan.cyhd.CN5.CNCHP.txt.bz2", package = "rCGH") cgh <- readAffyCytoScan(filePath, sampleName = "AffyScHD") cgh <- adjustSignal(cgh, nCores=1) cgh <- segmentCGH(cgh, nCores=1) cgh <- EMnormalize(cgh) ## Not run: # Interactive visalizations view(cgh) ## End(Not run)