plot_compare_profiles: Compare two 96 mutation profiles

Description Usage Arguments Value See Also Examples

View source: R/plot_compare_profiles.R

Description

Plots two 96 mutation profiles and their difference, reports the residual sum of squares (RSS).

Usage

1
2
3
4
5
6
7
8
9
plot_compare_profiles(
  profile1,
  profile2,
  profile_names = c("profile 1", "profile 2"),
  profile_ymax = 0.2,
  diff_ylim = c(-0.02, 0.02),
  colors = NA,
  condensed = FALSE
)

Arguments

profile1

First 96 mutation profile

profile2

Second 96 mutation profile

profile_names

Character vector with names of the mutations profiles used for plotting, default = c("profile 1", "profile 2")

profile_ymax

Maximum value of y-axis (relative contribution) for profile plotting. This can only be used to increase the y axis. If bars fall outside this limit, the maximum value is automatically increased. default = 0.2.

diff_ylim

Y-axis limits for profile difference plot, default = c(-0.02, 0.02)

colors

6 value color vector

condensed

More condensed plotting format. Default = F.

Value

96 spectrum plot of profile 1, profile 2 and their difference

See Also

mut_matrix, extract_signatures, plot_compare_indels, plot_compare_dbs, plot_compare_mbs

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
## See the 'mut_matrix()' example for how we obtained the following
## mutation matrix.
mut_mat <- readRDS(system.file("states/mut_mat_data.rds",
  package = "MutationalPatterns"
))

## Extracting signatures can be computationally intensive, so
## we use pre-computed data generated with the following command:
# nmf_res <- extract_signatures(mut_mat, rank = 2)

nmf_res <- readRDS(system.file("states/nmf_res_data.rds",
  package = "MutationalPatterns"
))

## Compare the reconstructed 96-profile of sample 1 with the original profile
## The same thing could be done with a reconstructed profile from signature refitting.
plot_compare_profiles(mut_mat[, 1],
  nmf_res$reconstructed[, 1],
  profile_names = c("Original", "Reconstructed")
)

## You could also compare regular mutation profiles with eachother.
plot_compare_profiles(
  mut_mat[, 1],
  mut_mat[, 2]
)


## You can also change the y limits.
## This can be done separately for the profiles and the different facets.
plot_compare_profiles(mut_mat[, 1],
  mut_mat[, 2],
  profile_ymax = 0.3,
  diff_ylim = c(-0.03, 0.03)
)

MutationalPatterns documentation built on Nov. 14, 2020, 2:03 a.m.