plotHeatmap | R Documentation |
Plots a heatmap of the copy number data. Each row is a cell and colums represent genomic positions.
plotHeatmap(
scCNA,
assay = "segment_ratios",
order_cells = NULL,
label = NULL,
label_colors = NULL,
group = NULL,
consensus = FALSE,
rounding_error = FALSE,
genes = NULL,
col = NULL,
row_split = NULL,
use_raster = TRUE,
raster_quality = 2,
n_threads = 1
)
scCNA |
The CopyKit object. |
assay |
String with the assay to pull data from to plot heatmap. |
order_cells |
A string with the desired method to order the cells within |
label |
A vector with the string names of the columns from
|
label_colors |
A named list with colors for the label annotation. Must match label length and have the same names as label |
group |
with the names of the columns from
|
consensus |
A boolean indicating if the consensus heatmap should be plotted. |
rounding_error |
A boolean indicating if the rounding error matrix should be plotted. |
genes |
A character vector with the HUGO symbol for genes to annotate on the heatmap. |
col |
A colorRamp2 vector that controls the color scale of the heatmap.
See |
row_split |
A string with the names of the columns from
|
use_raster |
Whether render the heatmap body as a raster image.
It helps to reduce file size when the matrix is huge. If number of rows or
columns is more than 2000, it is by default turned on. Note if cell_fun is
set, use_raster is enforced to be FALSE.
see |
raster_quality |
A value larger than 1. Larger values increase the file size. |
n_threads |
Number of threads passed on to |
order_cells: If order_cells argument is set to 'consensus_tree'
plotHeatmap
checks for the existence of a consensus matrix.
From the consensus matrix, a minimum evolution tree is built and cells are
ordered following the order of their respective groups from the tree.
If order_cells is 'hclust' cells are ordered according to hierarchical
clustering. 'hclust' calculation can be sped up by changing the parameter
'n_threads' if you have more threads available to use. If order_cells
is NULL the order of cells will be the same as the current order inside
the CopyKit object (colnames(CopyKit)).
label: A vector with the string names of the columns from
colData
for heatmap annotation. The
'label' argument can take as many columns as desired as argument as long
as they are elements from colData
.
label_colors: A named list, list element names must match column
names for colData
and list elements
must match the number of items present in the columns provided in argument
'label'. For example: to set colors for column 'outlier' containing
elements 'TRUE' or 'FALSE' a valid input would be:
'list(outlier = c('FALSE' = 'green', 'TRUE' = 'red))'.
Default colors are provided for 'superclones', 'subclones',
'is_aneuploid', and 'outlier' that can be override with 'label_colors'.
rounding_error: Must be used with assay = 'integer'.
plotHeatmap
will access the ploidies into colData(scCNA)$ploidy
that are generated from calcInteger
and scale rounded integer
values to the segment means. Later this scaled matrix will be subtracted
from the 'integer' assay from calcInteger
and the resulting
matrix from this subtraction will be plotted. Useful to visualize regions
of high rounding error. Such regions can indicate issues with the ploidy
scaling in use.
consensus: If set to TRUE, plotHeatmap
will search for
the consensus matrix in the slot consensus
and plot the
resulting matrix. Labels annotations can be added with the argument 'label'.
A ComplexHeatmap
object with a heatmap of copy number data
where the columns are the genomic positions and each row is a cell.
Darlan Conterno Minussi
Zuguang Gu, Roland Eils, Matthias Schlesner, Complex heatmaps reveal patterns and correlations in multidimensional genomic data, Bioinformatics, Volume 32, Issue 18, 15 September 2016, Pages 2847–2849, https://doi.org/10.1093/bioinformatics/btw313
calcInteger
.
copykit_obj <- copykit_example_filtered()
set.seed(1000)
copykit_obj <- copykit_obj[, sample(200)]
copykit_obj <- findClusters(copykit_obj)
copykit_obj <- calcConsensus(copykit_obj)
copykit_obj <- runConsensusPhylo(copykit_obj)
colData(copykit_obj)$section <- stringr::str_extract(
colData(copykit_obj)$sample,
"(L[0-9]+L[0-9]+|L[0-9]+)"
)
plotHeatmap(copykit_obj, label = c("section", "subclones"))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.