ContrastsPlotter: plot contrasts

ContrastsPlotterR Documentation

plot contrasts

Description

plot contrasts

plot contrasts

Public fields

contrastDF

data frame with contrasts

modelName

of column with model name

subject_Id

hierarchy key columns

prefix

default Contrasts - used to generate file names

diff

column with fold change differences

contrast

column with contrasts names, default "contrast"

volcano_spec

volcano plot specification

score_spec

score plot specification

histogram_spec

plot specification

fcthresh

fold change threshold

avg.abundance

name of column containing avg abundance values.

protein_annot

protein annotation

Methods

Public methods


Method new()

create Crontrast_Plotter

Usage
ContrastsPlotter$new(
  contrastDF,
  subject_Id,
  volcano = list(list(score = "FDR", thresh = 0.1)),
  histogram = list(list(score = "p.value", xlim = c(0, 1, 0.05)), list(score = "FDR",
    xlim = c(0, 1, 0.05))),
  score = list(list(score = "statistic", thresh = NULL)),
  fcthresh = 1,
  modelName = "modelName",
  diff = "diff",
  contrast = "contrast",
  avg.abundance = "avgAbd",
  protein_annot = NULL
)
Arguments
contrastDF

frame with contrast data

subject_Id

columns containing subject Identifier

volcano

which score to plot and which ablines to add.

histogram

which scores to plot and which range (x) should be shown.

score

score parameters

fcthresh

default 1 (log2 FC threshold)

modelName

name of column with model names

diff

fold change (difference) diff column

contrast

contrast column

avg.abundance

name of column with average abundance

protein_annot

add protein annotation (optional)


Method histogram()

plot histogram of selected scores (e.g. p-value, FDR, t-statistics)

Usage
ContrastsPlotter$histogram()

Method histogram_estimate()

plot histogram of effect size - difference between groups

Usage
ContrastsPlotter$histogram_estimate(binwidth = 0.05)
Arguments
binwidth

with of bin in histogram


Method histogram_diff()

plot histogram of differences (diff) fold change

Usage
ContrastsPlotter$histogram_diff(binwidth = 0.05)
Arguments
binwidth

with of bin in histogram


Method volcano()

volcano plots (fold change vs FDR)

Usage
ContrastsPlotter$volcano(
  colour,
  legend = TRUE,
  scales = c("fixed", "free", "free_x", "free_y"),
  min_score = 1e-04
)
Arguments
colour

column name with color information default modelName

legend

default TRUE

scales

default fixed facet_wrap, scales argument

min_score

replace p.values or FDR's smaller then min_score with min_score (default 0.0001).


Method volcano_plotly()

plotly volcano plots

Usage
ContrastsPlotter$volcano_plotly(
  colour,
  legend = TRUE,
  scales = c("fixed", "free", "free_x", "free_y"),
  min_score = 1e-04
)
Arguments
colour

column in contrast matrix with colour coding

legend

default TRUE

scales

default fixed facet_wrap, scales argument

min_score

replace p.values or FDR's smaller then min_score with min_score (default 0.0001).

Returns

list of ggplots


Method ma_plot()

ma plot

MA plot displays the effect size estimate as a function of the mean protein intensity across groups. Each dot represents an observed protein. Red horizontal lines represent the fold-change threshold.

Sometimes measured effects sizes (differences between samples groups) are biased by the signal intensity (here protein abundance). Such systematic effects can be explored using MA-plots.

Usage
ContrastsPlotter$ma_plot(fc, colour, legend = TRUE, rank = TRUE)
Arguments
fc

fold change abline

colour

column in contrast matrix with colour coding

legend

enable legend default TRUE

rank

default FALSE, if TRUE then rank of avgAbd is used.

Returns

ggplot


Method ma_plotly()

ma plotly

Usage
ContrastsPlotter$ma_plotly(fc, colour, legend = TRUE, rank = FALSE)
Arguments
fc

horizontal lines

colour

column in contrast matrix with colour coding.

legend

enable legend default TRUE

rank

default FALSE, if TRUE then rank of avgAbd is used.

Returns

list of ggplots


Method score_plot()

plot a score against the log2 fc e.g. t-statistic

Usage
ContrastsPlotter$score_plot(scorespec, colour, legend = TRUE)
Arguments
scorespec

list(score="statistics", fcthres = 2, thresh = 5)

colour

column with colour coding

legend

enable legend default TRUE

Returns

list of ggplots


Method score_plotly()

plot a score against the log2 fc e.g. t-statistic

Usage
ContrastsPlotter$score_plotly(scorespec, colour, legend = TRUE)
Arguments
scorespec

list(score="statistics", fcthres = 2, thresh = 5)

colour

column with colour coding

legend

enable legend default TRUE

Returns

list of ggplots


Method barplot_threshold()

shows the number of significant proteins per contrasts

Usage
ContrastsPlotter$barplot_threshold()
Returns

list which contains ggplots and summary tables


Method clone()

The objects of this class are cloneable with this method.

Usage
ContrastsPlotter$clone(deep = FALSE)
Arguments
deep

Whether to make a deep clone.

See Also

Other modelling: Contrasts, ContrastsMissing, ContrastsModerated, ContrastsProDA, ContrastsROPECA, ContrastsTable, INTERNAL_FUNCTIONS_BY_FAMILY, LR_test(), Model, build_model(), contrasts_fisher_exact(), get_anova_df(), get_complete_model_fit(), get_p_values_pbeta(), isSingular_lm(), linfct_all_possible_contrasts(), linfct_factors_contrasts(), linfct_from_model(), linfct_matrix_contrasts(), merge_contrasts_results(), model_analyse(), model_summary(), moderated_p_limma(), moderated_p_limma_long(), my_contest(), my_contrast(), my_contrast_V1(), my_contrast_V2(), my_glht(), pivot_model_contrasts_2_Wide(), plot_lmer_peptide_predictions(), sim_build_models_lm(), sim_build_models_lmer(), sim_make_model_lm(), sim_make_model_lmer(), strategy_lmer(), summary_ROPECA_median_p.scaled()

Other plotting: INTERNAL_FUNCTIONS_BY_FAMILY, UpSet_interaction_missing_stats(), UpSet_missing_stats(), medpolish_estimate_df(), missigness_histogram(), missingness_per_condition(), missingness_per_condition_cumsum(), plot_NA_heatmap(), plot_estimate(), plot_heatmap(), plot_heatmap_cor(), plot_hierarchies_add_quantline(), plot_hierarchies_boxplot_df(), plot_hierarchies_line(), plot_hierarchies_line_df(), plot_intensity_distribution_violin(), plot_pca(), plot_raster(), plot_sample_correlation(), plot_screeplot()

Examples


istar <- sim_lfq_data_peptide_config(Nprot = 100)
modelName <- "Model"
modelFunction <-
  strategy_lmer("abundance  ~ group_ + (1 | peptide_Id) + (1|sample)",
   model_name = modelName)
pepIntensity <- istar$data
config <- istar$config

mod <- build_model(
  pepIntensity,
  modelFunction,
  modelName = modelName,
  subject_Id = config$table$hierarchy_keys_depth())

 Contr <- c("group_A_vs_Ctrl" = "group_A - group_Ctrl",
 "group_B_vs_Ctrl" = "group_B - group_Ctrl"
 )
contrast <- prolfqua::Contrasts$new(mod,
  Contr)
contr <- contrast$get_contrasts()
cp <- ContrastsPlotter$new(contr,
  contrast$subject_Id,
  volcano = list(list(score = "FDR", thresh = 0.1)),
  histogram = list(list(score = "p.value", xlim = c(0,1,0.05)),
                 list(score = "FDR", xlim = c(0,1,0.05))),
  score =list(list(score = "statistic",  thresh = 5))
  )

stopifnot("plotly" %in% class(cp$volcano_plotly()$FDR))
stopifnot("ggplot" %in% class(cp$score_plot(legend=FALSE)$statistic))
p <- cp$histogram()
stopifnot("ggplot" %in% class(p$FDR))
stopifnot("ggplot" %in% class(p$p.value))
p <- cp$histogram_estimate()
stopifnot("ggplot" %in% class(p))
res <- cp$volcano()
stopifnot("ggplot" %in% class(res$FDR))
respltly <- cp$volcano_plotly()
stopifnot("plotly" %in% class(respltly$FDR))
stopifnot("ggplot" %in%  class(cp$ma_plot()))
stopifnot("plotly" %in%   class(cp$ma_plotly(rank=TRUE)))
res  <- cp$barplot_threshold()
stopifnot("ggplot" %in% class(res$FDR$plot))
stopifnot("ggplot" %in%  class(cp$histogram_diff()))


wolski/prolfqua documentation built on Dec. 4, 2024, 11:18 p.m.