Title: | Co-expression Modules identification Tool |
---|---|
Description: | The CEMiTool package unifies the discovery and the analysis of coexpression gene modules in a fully automatic manner, while providing a user-friendly html report with high quality graphs. Our tool evaluates if modules contain genes that are over-represented by specific pathways or that are altered in a specific sample group. Additionally, CEMiTool is able to integrate transcriptomic data with interactome information, identifying the potential hubs on each network. |
Authors: | Pedro Russo [aut], Gustavo Ferreira [aut], Matheus Bürger [aut], Lucas Cardozo [aut], Diogenes Lima [aut], Thiago Hirata [aut], Melissa Lever [aut], Helder Nakaya [aut, cre] |
Maintainer: | Helder Nakaya <[email protected]> |
License: | GPL-3 |
Version: | 1.31.0 |
Built: | 2024-10-30 04:35:58 UTC |
Source: | https://github.com/bioc/CEMiTool |
This function takes a CEMiTool
object containing expression values
and returns a CEMiTool object with an adjacency matrix in the adjacency slot.
adj_data(cem, ...) ## S4 method for signature 'CEMiTool' adj_data(cem) adj_data(cem) <- value ## S4 replacement method for signature 'CEMiTool' adj_data(cem) <- value
adj_data(cem, ...) ## S4 method for signature 'CEMiTool' adj_data(cem) adj_data(cem) <- value ## S4 replacement method for signature 'CEMiTool' adj_data(cem) <- value
cem |
Object of class |
... |
Optional parameters. |
value |
Object of class |
Object of class matrix
with adjacency values or object of class CEMiTool
.
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Calculate adjacency matrix with example beta value 8 cem <- get_adj(cem, beta=8) # Return adjacency matrix adj <- adj_data(cem) # Check result adj[1:5, 1:5] # Set adjacency matrix to CEMiTool object adj_data(cem) <- adj
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Calculate adjacency matrix with example beta value 8 cem <- get_adj(cem, beta=8) # Return adjacency matrix adj <- adj_data(cem) # Check result adj[1:5, 1:5] # Set adjacency matrix to CEMiTool object adj_data(cem) <- adj
This object can be used as input for CEMiTool functions. Data used are from
expr
and sample_annot
.
data(cem)
data(cem)
An object of class CEMiTool
# Get example CEMiTool object data(cem) cem
# Get example CEMiTool object data(cem) cem
Defines co-expression modules and runs several different analyses.
cemitool( expr, annot, gmt, interactions, filter = TRUE, filter_pval = 0.1, apply_vst = FALSE, n_genes, eps = 0.1, cor_method = c("pearson", "spearman"), cor_function = "cor", network_type = "unsigned", tom_type = "signed", set_beta = NULL, force_beta = FALSE, sample_name_column = "SampleName", class_column = "Class", merge_similar = TRUE, rank_method = "mean", ora_pval = 0.05, gsea_scale = TRUE, gsea_min_size = 15, gsea_max_size = 1000, min_ngen = 30, diss_thresh = 0.8, plot = TRUE, plot_diagnostics = TRUE, order_by_class = TRUE, center_func = "mean", directed = FALSE, verbose = FALSE )
cemitool( expr, annot, gmt, interactions, filter = TRUE, filter_pval = 0.1, apply_vst = FALSE, n_genes, eps = 0.1, cor_method = c("pearson", "spearman"), cor_function = "cor", network_type = "unsigned", tom_type = "signed", set_beta = NULL, force_beta = FALSE, sample_name_column = "SampleName", class_column = "Class", merge_similar = TRUE, rank_method = "mean", ora_pval = 0.05, gsea_scale = TRUE, gsea_min_size = 15, gsea_max_size = 1000, min_ngen = 30, diss_thresh = 0.8, plot = TRUE, plot_diagnostics = TRUE, order_by_class = TRUE, center_func = "mean", directed = FALSE, verbose = FALSE )
expr |
Gene expression |
annot |
Sample annotation |
gmt |
A data.frame containing two columns, one with pathways and one with genes |
interactions |
A data.frame containing two columns with gene names. |
filter |
Logical. If TRUE, will filter expression data. |
filter_pval |
P-value threshold for filtering. Default |
apply_vst |
Logical. If TRUE, will apply Variance Stabilizing Transform before filtering genes.
Currently ignored if parameter |
n_genes |
Number of genes left after filtering. |
eps |
A value for accepted R-squared interval between subsequent beta values. Default is 0.1. |
cor_method |
A character string indicating which correlation coefficient is
to be computed. One of |
cor_function |
A character string indicating the correlation function to be used. Supported functions are
currently 'cor' and 'bicor'. Default is |
network_type |
A character string indicating if network type should be computed
as |
tom_type |
A character string indicating if the TOM type should be computed
as |
set_beta |
A value to override the automatically selected beta value. Default is NULL. |
force_beta |
Whether or not to automatically force a beta value based on number of samples. Default is FALSE. |
sample_name_column |
A character string indicating the sample column name of the annotation table. |
class_column |
A character string indicating the class column name of the annotation table. |
merge_similar |
Logical. If |
rank_method |
Character string indicating how to rank genes. Either "mean" (the default) or "median". |
ora_pval |
P-value for overrepresentation analysis. Default |
gsea_scale |
If TRUE, apply z-score transformation for GSEA analysis. Default is |
gsea_min_size |
Minimum size of gene sets for GSEA analysis. Default is |
gsea_max_size |
Maximum size of gene sets for GSEA analysis. Default is |
min_ngen |
Minimal number of genes per submodule. Default |
diss_thresh |
Module merging correlation threshold for eigengene similarity.
Default |
plot |
Logical. If |
plot_diagnostics |
Logical. If |
order_by_class |
Logical. If |
center_func |
Character string indicating the centrality measure to show in the plot. Either 'mean' (the default) or 'median'. |
directed |
Logical. If |
verbose |
Logical. If |
Object of class CEMiTool
# Get example expression data data(expr0) # Run CEMiTool analyses cem <- cemitool(expr=expr0) # Run CEMiTool applying Variance Stabilizing Transformation to data cem <- cemitool(expr=expr0, apply_vst=TRUE) # Run CEMiTool with additional processing messages cem <- cemitool(expr=expr0, verbose=TRUE) ## Not run: # Run full CEMiTool analysis ## Get example sample annotation data data(sample_annot) ## Read example pathways file gmt_fname <- system.file("extdata", "pathways.gmt", package = "CEMiTool") gmt_in <- read_gmt(gmt_fname) ## Get example interactions file int_df <- read.delim(system.file("extdata", "interactions.tsv", package = "CEMiTool")) ## Run CEMiTool cem <- cemitool(expr=expr0, annot=sample_annot, gmt=gmt_in, interactions=int_df, verbose=TRUE, plot=TRUE) # Create report as html file generate_report(cem, directory = "./Report", output_format="html_document") # Write analysis results into files write_files(cem, directory="./Tables", force=TRUE) # Save all plots save_plots(cem, "all", directory="./Plots") ## End(Not run)
# Get example expression data data(expr0) # Run CEMiTool analyses cem <- cemitool(expr=expr0) # Run CEMiTool applying Variance Stabilizing Transformation to data cem <- cemitool(expr=expr0, apply_vst=TRUE) # Run CEMiTool with additional processing messages cem <- cemitool(expr=expr0, verbose=TRUE) ## Not run: # Run full CEMiTool analysis ## Get example sample annotation data data(sample_annot) ## Read example pathways file gmt_fname <- system.file("extdata", "pathways.gmt", package = "CEMiTool") gmt_in <- read_gmt(gmt_fname) ## Get example interactions file int_df <- read.delim(system.file("extdata", "interactions.tsv", package = "CEMiTool")) ## Run CEMiTool cem <- cemitool(expr=expr0, annot=sample_annot, gmt=gmt_in, interactions=int_df, verbose=TRUE, plot=TRUE) # Create report as html file generate_report(cem, directory = "./Report", output_format="html_document") # Write analysis results into files write_files(cem, directory="./Tables", force=TRUE) # Save all plots save_plots(cem, "all", directory="./Plots") ## End(Not run)
An S4 class to represent the CEMiTool analysis.
expression
Gene expression data.frame
.
sample_annotation
Sample annotation data.frame
.
fit_indices
data.frame
containing scale-free model fit,
soft-threshold and network parameters.
selected_genes
Character vector
containing the names of genes
selected for analysis
module
Genes in modules information data.frame
.
enrichment
list
with modules enrichment results for sample classes.
ora
Over-representation analysis results data.frame
.
interactions
list
containing gene interactions present in
modules.
interaction_plot
list of ggplot graphs with module gene interactions.
profile_plot
list of ggplot graphs with gene expression profile per module.
enrichment_plot
ggplot graph for enrichment analysis results.
beta_r2_plot
ggplot graph with scale-free topology fit results for each soft-threshold.
mean_k_plot
ggplot graph with mean network connectivity.
barplot_ora
list of ggplot graphs with over-representation analysis results per module.
sample_tree_plot
gtable containing sample dendrogram with class labels and clinical data (if available in sample_annotation(cem)).
mean_var_plot
Mean x variance scatterplot.
hist_plot
Expression histogram.
qq_plot
Quantile-quantile plot.
sample_name_column
character string containing the name of the column with sample names in the annotation file.
class_column
character string containing the name of the column with class names in the annotation file.
mod_colors
character vector
containing colors associated with each network module.
parameters
list
containing analysis parameters.
adjacency
matrix
containing gene adjacency values based on correlation
# Get example expression data data(expr0) # Initialize CEMiTool object with expression cem <- new("CEMiTool", expression=expr0)
# Get example expression data data(expr0) # Initialize CEMiTool object with expression cem <- new("CEMiTool", expression=expr0)
Creates report for identifying potential problems with data.
diagnostic_report(cem, ...) ## S4 method for signature 'CEMiTool' diagnostic_report( cem, title = "Diagnostics", directory = "./Reports/Diagnostics", force = FALSE, ... )
diagnostic_report(cem, ...) ## S4 method for signature 'CEMiTool' diagnostic_report( cem, title = "Diagnostics", directory = "./Reports/Diagnostics", force = FALSE, ... )
cem |
Object of class |
... |
parameters to rmarkdown::render |
title |
Character string with the title of the report. |
directory |
Directory name for results. |
force |
If the directory exists, execution will not stop. |
An HTML file with an interactive diagnostic report.
Retrieve and set expression attribute
expr_data(cem, ...) ## S4 method for signature 'CEMiTool' expr_data(cem, filter = TRUE, apply_vst = FALSE, filter_pval = 0.1, ...) expr_data(cem) <- value ## S4 replacement method for signature 'CEMiTool' expr_data(cem) <- value
expr_data(cem, ...) ## S4 method for signature 'CEMiTool' expr_data(cem, filter = TRUE, apply_vst = FALSE, filter_pval = 0.1, ...) expr_data(cem) <- value ## S4 replacement method for signature 'CEMiTool' expr_data(cem) <- value
cem |
Object of class |
... |
Additional parameters to |
filter |
logical. If TRUE, retrieves filtered expression data (Default: TRUE) |
apply_vst |
logical. If TRUE, applies variance stabilizing transformation to expression data (Default: FALSE) |
filter_pval |
logical. Threshold for filter p-value. Ignored if filter = FALSE (Default: 0.1) |
value |
Object of class |
Object of class data.frame
with gene expression data
# Initialize an empty CEMiTool object cem <- new_cem() # Get example expression data data(expr0) # Add expression file to CEMiTool object expr_data(cem) <- expr0 # Check expression file head(expr_data(cem))
# Initialize an empty CEMiTool object cem <- new_cem() # Get example expression data data(expr0) # Add expression file to CEMiTool object expr_data(cem) <- expr0 # Check expression file head(expr_data(cem))
Modified data from a yellow fever vaccination study by Querec et al, 2009. In order to reduce package size, only the 4000 genes with the highest variance were selected for this dataset.
data(expr0)
data(expr0)
An object of class data.frame
Querec TD, Akondy RS, Lee EK, Cao W et al. Systems biology approach predicts immunogenicity of the yellow fever vaccine in humans. Nat Immunol 2009 Jan;10(1):116-25. PMID: 19029902 PubMed
data(expr0) # Run CEMiTool analysis ## Not run: cemitool(expr0)
data(expr0) # Run CEMiTool analysis ## Not run: cemitool(expr0)
Filter a gene expression data.frame
filter_genes(expr, pct = 0.75, apply_vst = FALSE)
filter_genes(expr, pct = 0.75, apply_vst = FALSE)
expr |
A data.frame containing expression data |
pct |
Percentage of most expressed genes to keep. |
apply_vst |
Logical. If TRUE, will apply variance stabilizing transform before filtering data. |
This function takes a gene expression data.frame and applies a
percentage filter (keeps the pct
) most expressed genes). If
apply_vst
is TRUE, a variance stabilizing transformation is applied on
gene expression values as long as mean and variance values have a Spearman's
rho of over 0.5. This transformation is intended to remove this dependence
between the parameters. One should then apply the select_genes
function to get significant genes.
A data.frame containing filtered expression data
# Get example expression data data(expr0) # Filter genes expr_f <- filter_genes(expr0) # Check selected genes expr_f[1:5, 1:5] # Filter genes and apply variance stabilizing transformation expr_f2 <- filter_genes(expr0, apply_vst=TRUE) # Check results expr_f2[1:5, 1:5] # Selected genes selected <- select_genes(expr_f2) # Get data.frame with only selected genes expr_s <- expr_f2[selected, ] # Check results expr_s[1:5, 1:5]
# Get example expression data data(expr0) # Filter genes expr_f <- filter_genes(expr0) # Check selected genes expr_f[1:5, 1:5] # Filter genes and apply variance stabilizing transformation expr_f2 <- filter_genes(expr0, apply_vst=TRUE) # Check results expr_f2[1:5, 1:5] # Selected genes selected <- select_genes(expr_f2) # Get data.frame with only selected genes expr_s <- expr_f2[selected, ] # Check results expr_s[1:5, 1:5]
Defines co-expression modules
find_modules(cem, ...) ## S4 method for signature 'CEMiTool' find_modules( cem, cor_method = c("pearson", "spearman"), cor_function = "cor", eps = 0.1, set_beta = NULL, force_beta = FALSE, min_ngen = 20, merge_similar = TRUE, network_type = "unsigned", tom_type = "signed", diss_thresh = 0.8, verbose = FALSE )
find_modules(cem, ...) ## S4 method for signature 'CEMiTool' find_modules( cem, cor_method = c("pearson", "spearman"), cor_function = "cor", eps = 0.1, set_beta = NULL, force_beta = FALSE, min_ngen = 20, merge_similar = TRUE, network_type = "unsigned", tom_type = "signed", diss_thresh = 0.8, verbose = FALSE )
cem |
Object of class |
... |
Optional parameters. |
cor_method |
A character string indicating which correlation coefficient
is to be computed. Default |
cor_function |
A character string indicating the correlation function to be used. Default |
eps |
A value for accepted R-squared interval between subsequent beta values. Default is 0.1. |
set_beta |
A value to override the automatically selected beta value. Default is NULL. |
force_beta |
Whether or not to automatically force a beta value based on number of samples. Default is FALSE. |
min_ngen |
Minimal number of genes per submodule. Default |
merge_similar |
Logical. If |
network_type |
A character string indicating to use either "unsigned" (default) or "signed" networks. Default |
tom_type |
A character string indicating to use either "unsigned" or "signed" (default) TOM similarity measure. |
diss_thresh |
Module merging correlation threshold for eigengene similarity.
Default |
verbose |
Logical. If |
Object of class CEMiTool
# Get example expression data data(expr0) # Initialize CEMiTool object with expression cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Define network modules cem <- find_modules(cem) # Check results head(module_genes(cem))
# Get example expression data data(expr0) # Initialize CEMiTool object with expression cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Define network modules cem <- find_modules(cem) # Check results head(module_genes(cem))
Retrieve scale-free model fit data
fit_data(cem) ## S4 method for signature 'CEMiTool' fit_data(cem)
fit_data(cem) ## S4 method for signature 'CEMiTool' fit_data(cem)
cem |
Object of class |
Object of class data.frame
# Get example CEMiTool object data(cem) # Get modules and beta data cem <- find_modules(cem) # Get fit data fit_data(cem)
# Get example CEMiTool object data(cem) # Get modules and beta data cem <- find_modules(cem) # Get fit data fit_data(cem)
Creates report for CEMiTool results
generate_report(cem, ...) ## S4 method for signature 'CEMiTool' generate_report( cem, max_rows_ora = 50, title = "Report", directory = "./Reports/Report", force = FALSE, ... )
generate_report(cem, ...) ## S4 method for signature 'CEMiTool' generate_report( cem, max_rows_ora = 50, title = "Report", directory = "./Reports/Report", force = FALSE, ... )
cem |
Object of class |
... |
parameters to rmarkdown::render |
max_rows_ora |
maximum number of rows in Over Representation Analysis table results |
title |
Character string with the title of the report. |
directory |
Directory name for results. |
force |
If the directory exists, execution will not stop. |
An HTML file with an interactive report of CEMiTool analyses.
## Not run: # Get example CEMiTool object data(cem) generate_report(cem, output_format=c("pdf_document", "html_document")) ## End(Not run)
## Not run: # Get example CEMiTool object data(cem) generate_report(cem, output_format=c("pdf_document", "html_document")) ## End(Not run)
This function takes a CEMiTool
object
and returns an adjacency matrix.
get_adj(cem, ...) ## S4 method for signature 'CEMiTool' get_adj( cem, beta, network_type = "unsigned", cor_function = "cor", cor_method = "pearson" )
get_adj(cem, ...) ## S4 method for signature 'CEMiTool' get_adj( cem, beta, network_type = "unsigned", cor_function = "cor", cor_method = "pearson" )
cem |
Object of class |
... |
Optional parameters. |
beta |
Selected soft-threshold value |
network_type |
A character string indicating to use either "unsigned"
(default) or "signed" networks. Default |
cor_function |
A character string indicating the correlation function
to be used. Default |
cor_method |
A character string indicating which correlation
coefficient is to be computed. Default |
Object of class CEMiTool
with adjacency data
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Calculate adjacency matrix with example beta value 8 cem <- get_adj(cem, beta=8) # Check results adj <- adj_data(cem) adj[1:5, 1:5]
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Calculate adjacency matrix with example beta value 8 cem <- get_adj(cem, beta=8) # Check results adj <- adj_data(cem) adj[1:5, 1:5]
This function takes the input parameters from find_modules and calculates the WGCNA soft-threshold parameters and returns them.
get_beta_data(cem, ...) ## S4 method for signature 'CEMiTool' get_beta_data( cem, network_type = "unsigned", cor_function = "cor", cor_method = "pearson", verbose = FALSE )
get_beta_data(cem, ...) ## S4 method for signature 'CEMiTool' get_beta_data( cem, network_type = "unsigned", cor_function = "cor", cor_method = "pearson", verbose = FALSE )
cem |
A CEMiTool object containing expression data |
... |
Optional parameters. |
network_type |
A character string indicating to use either "unsigned"
(default) or "signed" networks. Default |
cor_function |
A character string indicating the correlation function
to be used. Default |
cor_method |
A character string indicating which correlation
coefficient is to be computed. Default |
verbose |
Logical. If |
A list containing the soft-threshold selected by WGCNA and scale-free model parameters
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Get beta data beta_data <- get_beta_data(cem)
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Get beta data beta_data <- get_beta_data(cem)
This function takes a CEMiTool object with beta data and returns a vector containing the chosen beta and corresponding R squared value.
get_cemitool_r2_beta(cem, ...) ## S4 method for signature 'CEMiTool' get_cemitool_r2_beta(cem, eps = 0.1)
get_cemitool_r2_beta(cem, ...) ## S4 method for signature 'CEMiTool' get_cemitool_r2_beta(cem, eps = 0.1)
cem |
A CEMiTool object containing the fit_indices slot |
... |
Optional parameters. |
eps |
A value indicating the accepted interval between successive values of R squared to use to calculate the selected beta. Default: 0.1. |
A vector containing R squared value and the chosen beta parameter.
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Get modules and beta data cem <- find_modules(cem) # Get CEMiTool R2 and beta values get_cemitool_r2_beta(cem)
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Get modules and beta data cem <- find_modules(cem) # Get CEMiTool R2 and beta values get_cemitool_r2_beta(cem)
This function takes a CEMiTool object and returns the network connectivity.
get_connectivity(cem, ...) ## S4 method for signature 'CEMiTool' get_connectivity(cem, beta)
get_connectivity(cem, ...) ## S4 method for signature 'CEMiTool' get_connectivity(cem, beta)
cem |
Object of class |
... |
Optional parameters. |
beta |
A soft-thresholding value to be used for the network. |
The value of the network's connectivity.
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Get modules and beta data cem <- find_modules(cem) # Get network connectivity with example beta value 8 get_connectivity(cem, beta=8)
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Get modules and beta data cem <- find_modules(cem) # Get network connectivity with example beta value 8 get_connectivity(cem, beta=8)
Returns n
genes in each module with high connectivity.
get_hubs(cem, ...) ## S4 method for signature 'CEMiTool' get_hubs(cem, n = 5, method = "adjacency")
get_hubs(cem, ...) ## S4 method for signature 'CEMiTool' get_hubs(cem, n = 5, method = "adjacency")
cem |
Object of class |
... |
Optional parameters. |
n |
Number of hubs to return from each module. If "all", returns all genes in decreasing order of connectivity. Default: 5. |
method |
Method for hub calculation. Either "adjacency" or "kME". Default: "adjacency" |
A list
containing hub genes for each module and the value of
the calculated method.
# Get example CEMiTool object data(cem) # Get module hubs hubs <- get_hubs(cem, n=10, "adjacency")
# Get example CEMiTool object data(cem) # Get module hubs hubs <- get_hubs(cem, n=10, "adjacency")
This function takes a CEMiTool object with expression and a vector of numeric labels to merge similar modules.
get_merged_mods(cem, ...) ## S4 method for signature 'CEMiTool' get_merged_mods(cem, mods, diss_thresh = 0.8, verbose = FALSE)
get_merged_mods(cem, ...) ## S4 method for signature 'CEMiTool' get_merged_mods(cem, mods, diss_thresh = 0.8, verbose = FALSE)
cem |
Object of class |
... |
Optional parameters. |
mods |
A vector containing numeric labels for each module gene |
diss_thresh |
A threshold of dissimilarity for modules. Default is 0.8. |
verbose |
Logical. If |
Numeric labels assigning genes to modules.
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Calculate adjacency matrix with example beta value 8 cem <- get_adj(cem, beta=8) # Get modules mods <- get_mods(cem) # Merge similar modules merged_mods <- get_merged_mods(cem, mods)
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Calculate adjacency matrix with example beta value 8 cem <- get_adj(cem, beta=8) # Get modules mods <- get_mods(cem) # Merge similar modules merged_mods <- get_merged_mods(cem, mods)
This function takes a CEMiTool
object containing an adjacency matrix
together with the given network parameters, and returns the given
co-expression modules.
get_mods(cem, ...) ## S4 method for signature 'CEMiTool' get_mods( cem, cor_function = "cor", cor_method = "pearson", tom_type = "signed", min_ngen = 20 )
get_mods(cem, ...) ## S4 method for signature 'CEMiTool' get_mods( cem, cor_function = "cor", cor_method = "pearson", tom_type = "signed", min_ngen = 20 )
cem |
Object of class |
... |
Optional parameters. |
cor_function |
A character string indicating the correlation function
to be used. Default |
cor_method |
A character string indicating which correlation
coefficient is to be computed. Default |
tom_type |
A character string indicating to use either "unsigned" or "signed" (default) TOM similarity measure. |
min_ngen |
Minimal number of genes per module (Default: 20). |
Numeric labels assigning genes to modules.
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Calculate adjacency matrix with example beta value 8 cem <- get_adj(cem, beta=8) # Get module labels mods <- get_mods(cem)
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Calculate adjacency matrix with example beta value 8 cem <- get_adj(cem, beta=8) # Get module labels mods <- get_mods(cem)
This function takes a CEMiTool object and returns the phi parameter.
get_phi(cem, ...) ## S4 method for signature 'CEMiTool' get_phi(cem)
get_phi(cem, ...) ## S4 method for signature 'CEMiTool' get_phi(cem)
cem |
A CEMiTool object containing the fit_indices slot |
... |
Optional parameters. |
The phi parameter
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Get modules and beta data cem <- find_modules(cem) # Get phi get_phi(cem)
# Get example expression data data(expr0) # Initialize new CEMiTool object with expression data cem <- new_cem(expr0, filter=TRUE, apply_vst=FALSE) # Get modules and beta data cem <- find_modules(cem) # Get phi get_phi(cem)
Retrieve Gene Set Enrichment Analysis (GSEA) results
gsea_data(cem) ## S4 method for signature 'CEMiTool' gsea_data(cem)
gsea_data(cem) ## S4 method for signature 'CEMiTool' gsea_data(cem)
cem |
Object of class |
Object of class list
with GSEA data
# Get example CEMiTool object data(cem) # Look at example annotation file sample_annotation(cem) # Run GSEA on network modules cem <- mod_gsea(cem) # Check results gsea_data(cem)
# Get example CEMiTool object data(cem) # Look at example annotation file sample_annotation(cem) # Run GSEA on network modules cem <- mod_gsea(cem) # Check results gsea_data(cem)
Retrieve and set interaction data to CEMiTool object
interactions_data(cem, ...) ## S4 method for signature 'CEMiTool' interactions_data(cem) interactions_data(cem) <- value ## S4 replacement method for signature 'CEMiTool' interactions_data(cem) <- value
interactions_data(cem, ...) ## S4 method for signature 'CEMiTool' interactions_data(cem) interactions_data(cem) <- value ## S4 replacement method for signature 'CEMiTool' interactions_data(cem) <- value
cem |
Object of class |
... |
parameters for igraph::graph_from_data_frame |
value |
a data.frame or matrix containing two columns |
Object of class CEMiTool
# Get example CEMiTool object data(cem) # Read example interactions data int_df <- read.delim(system.file("extdata", "interactions.tsv", package = "CEMiTool")) # Insert interactions data interactions_data(cem) <- int_df # Check interactions data interactions_data(cem)
# Get example CEMiTool object data(cem) # Read example interactions data int_df <- read.delim(system.file("extdata", "interactions.tsv", package = "CEMiTool")) # Insert interactions data interactions_data(cem) <- int_df # Check interactions data interactions_data(cem)
Retrieve and set mod_colors attribute
mod_colors(cem) ## S4 method for signature 'CEMiTool' mod_colors(cem) mod_colors(cem) <- value ## S4 replacement method for signature 'CEMiTool,character' mod_colors(cem) <- value
mod_colors(cem) ## S4 method for signature 'CEMiTool' mod_colors(cem) mod_colors(cem) <- value ## S4 replacement method for signature 'CEMiTool,character' mod_colors(cem) <- value
cem |
Object of class |
value |
a character vector containing colors for each module. Names should match with module names |
A vector with color names.
# Get example CEMiTool object data(cem) # See module colors mod_colors(cem)
# Get example CEMiTool object data(cem) # See module colors mod_colors(cem)
Get the number of genes in modules in a CEMiTool object
mod_gene_num(cem, module = NULL) ## S4 method for signature 'CEMiTool' mod_gene_num(cem, module = NULL)
mod_gene_num(cem, module = NULL) ## S4 method for signature 'CEMiTool' mod_gene_num(cem, module = NULL)
cem |
Object of class |
module |
Default is NULL. If a character string designating a module is given, the number of genes in that module is returned instead. |
The number of genes in module(s)
# Get example CEMiTool object data(cem) # Get the number of genes in modules mod_gene_num(cem) # Get the number of genes in module M1 mod_gene_num(cem, "M1")
# Get example CEMiTool object data(cem) # Get the number of genes in modules mod_gene_num(cem) # Get the number of genes in module M1 mod_gene_num(cem, "M1")
Perfoms Gene Set Enrichment Analysis (GSEA) for each co-expression module found.
mod_gsea(cem, ...) ## S4 method for signature 'CEMiTool' mod_gsea( cem, gsea_scale = TRUE, rank_method = "mean", gsea_min_size = 15, gsea_max_size = 1000, verbose = FALSE )
mod_gsea(cem, ...) ## S4 method for signature 'CEMiTool' mod_gsea( cem, gsea_scale = TRUE, rank_method = "mean", gsea_min_size = 15, gsea_max_size = 1000, verbose = FALSE )
cem |
Object of class |
... |
Optional parameters. |
gsea_scale |
If TRUE, transform data using z-score transformation. Default: TRUE |
rank_method |
Character string indicating how to rank genes. Either "mean" (the default) or "median". |
gsea_min_size |
Minimum gene set size (Default: 15). |
gsea_max_size |
Maximum gene set size (Default: 1000). |
verbose |
logical. Report analysis steps. |
GSEA results.
# Get example CEMiTool object data(cem) # Look at example annotation file sample_annotation(cem) # Run GSEA on network modules cem <- mod_gsea(cem) # Check results gsea_data(cem)
# Get example CEMiTool object data(cem) # Look at example annotation file sample_annotation(cem) # Run GSEA on network modules cem <- mod_gsea(cem) # Check results gsea_data(cem)
Get module names in a CEMiTool object
mod_names(cem, include_NC = TRUE) ## S4 method for signature 'CEMiTool' mod_names(cem, include_NC = TRUE)
mod_names(cem, include_NC = TRUE) ## S4 method for signature 'CEMiTool' mod_names(cem, include_NC = TRUE)
cem |
Object of class |
include_NC |
Logical. Whether or not to include "Not.Correlated"
module. Defaults to |
Module names
# Get example CEMiTool object data(cem) # Get module names mod_names(cem)
# Get example CEMiTool object data(cem) # Get module names mod_names(cem)
Performs overrepresentation analysis for each co-expression module found.
mod_ora(cem, ...) ## S4 method for signature 'CEMiTool' mod_ora(cem, gmt, verbose = FALSE)
mod_ora(cem, ...) ## S4 method for signature 'CEMiTool' mod_ora(cem, gmt, verbose = FALSE)
cem |
Object of class |
... |
Optional parameters. |
gmt |
Object of class |
verbose |
logical. Report analysis steps. |
Object of class CEMiTool
# Get example CEMiTool object data(cem) # Read gmt file gmt <- read_gmt(system.file('extdata', 'pathways.gmt', package='CEMiTool')) # Run module overrepresentation analysis cem <- mod_ora(cem, gmt) # Check results head(ora_data(cem))
# Get example CEMiTool object data(cem) # Read gmt file gmt <- read_gmt(system.file('extdata', 'pathways.gmt', package='CEMiTool')) # Run module overrepresentation analysis cem <- mod_ora(cem, gmt) # Check results head(ora_data(cem))
Summarizes modules using mean or eigengene expression.
mod_summary(cem, ...) ## S4 method for signature 'CEMiTool' mod_summary(cem, method = c("mean", "median", "eigengene"), verbose = FALSE)
mod_summary(cem, ...) ## S4 method for signature 'CEMiTool' mod_summary(cem, method = c("mean", "median", "eigengene"), verbose = FALSE)
cem |
Object of class |
... |
Optional parameters. |
method |
A character string indicating which summarization method is to be used. Can be 'eigengene', 'mean' or 'median'. Default is 'mean'. |
verbose |
Logical. If |
A data.frame
with summarized values.
# Get example CEMiTool object data(cem) # Summarize results mod_summary <- mod_summary(cem)
# Get example CEMiTool object data(cem) # Summarize results mod_summary <- mod_summary(cem)
Get the module genes in a CEMiTool object
module_genes(cem, module = NULL) ## S4 method for signature 'CEMiTool' module_genes(cem, module = NULL)
module_genes(cem, module = NULL) ## S4 method for signature 'CEMiTool' module_genes(cem, module = NULL)
cem |
Object of class |
module |
A character string with the name of the module of which
genes are to be returned. Defaults to |
Object of class data.frame
containing genes and their
respective module
# Get example CEMiTool object data(cem) # Get the module genes module_genes(cem) # Get genes for module M1 module_genes(cem, module="M1")
# Get example CEMiTool object data(cem) # Get the module genes module_genes(cem) # Get genes for module M1 module_genes(cem, module="M1")
Create a CEMiTool object
new_cem( expr = data.frame(), sample_annot = data.frame(), sample_name_column = "SampleName", class_column = "Class", filter = TRUE, apply_vst = FALSE, filter_pval = 0.1 )
new_cem( expr = data.frame(), sample_annot = data.frame(), sample_name_column = "SampleName", class_column = "Class", filter = TRUE, apply_vst = FALSE, filter_pval = 0.1 )
expr |
Object of class |
sample_annot |
Object of |
sample_name_column |
A string specifying the column to be used as sample identification. Default: "SampleName". |
class_column |
A string specifying the column to be used as a grouping factor for samples. Default: "Class" |
filter |
Logical. Used to define if posterior functions should use filtered expression data or not (Default: TRUE) |
apply_vst |
Logical. Used to define if posterior functions should use a
variance stabilizing transformation on expression data before analyses. Only
valid if argument |
filter_pval |
logical. Threshold for filter p-value. Ignored if filter = FALSE (Default: 0.1) |
Object of class CEMiTool
# Create new CEMiTool object cem <- new_cem() # Create new CEMiTool object with expression and sample_annotation data data(expr0) data(sample_annot) cem <- new_cem(expr0, sample_annot, "SampleName", "Class") # Equivalent to a call to new() cem2 <- new("CEMiTool", expression=expr0, sample_annotation=sample_annot) identical(cem, cem2)
# Create new CEMiTool object cem <- new_cem() # Create new CEMiTool object with expression and sample_annotation data data(expr0) data(sample_annot) cem <- new_cem(expr0, sample_annot, "SampleName", "Class") # Equivalent to a call to new() cem2 <- new("CEMiTool", expression=expr0, sample_annotation=sample_annot) identical(cem, cem2)
Get the number of modules in a CEMiTool object
nmodules(cem) ## S4 method for signature 'CEMiTool' nmodules(cem)
nmodules(cem) ## S4 method for signature 'CEMiTool' nmodules(cem)
cem |
Object of class |
number of modules
# Get example CEMiTool object data(cem) # Get the number of modules nmodules(cem)
# Get example CEMiTool object data(cem) # Get the number of modules nmodules(cem)
Retrieve over representation analysis (ORA) results
ora_data(cem) ## S4 method for signature 'CEMiTool' ora_data(cem)
ora_data(cem) ## S4 method for signature 'CEMiTool' ora_data(cem)
cem |
Object of class |
This function returns the results of the mod_ora
function on the
CEMiTool
object. The ID column corresponds to pathways in the gmt file for which
genes in the modules were enriched. The Count column shows the number of genes in the
module that are enriched for each pathway. The GeneRatio column shows the proportion of
genes in the module enriched for a given pathway out of all the genes in the module
enriched for any given pathway. The BgRatio column shows the proportion of genes in a
given pathway out of all the genes in the gmt file. For more details, please refer to
the clusterProfiler
package documentation.
Object of class data.frame
with ORA data
Guangchuang Yu, Li-Gen Wang, Yanyan Han, Qing-Yu He. clusterProfiler: an R package for comparing biological themes among gene clusters. OMICS: A Journal of Integrative Biology. 2012, 16(5):284-287.
# Get example CEMiTool object data(cem) # Read gmt file gmt <- read_gmt(system.file('extdata', 'pathways.gmt', package='CEMiTool')) # Run module overrepresentation analysis cem <- mod_ora(cem, gmt) # Check results head(ora_data(cem))
# Get example CEMiTool object data(cem) # Read gmt file gmt <- read_gmt(system.file('extdata', 'pathways.gmt', package='CEMiTool')) # Run module overrepresentation analysis cem <- mod_ora(cem, gmt) # Check results head(ora_data(cem))
Creates a graph showing each possible soft-threshold value and its corresponding R squared value
plot_beta_r2(cem, ...) ## S4 method for signature 'CEMiTool' plot_beta_r2(cem, plot_title = "Scale independence (beta selection)")
plot_beta_r2(cem, ...) ## S4 method for signature 'CEMiTool' plot_beta_r2(cem, plot_title = "Scale independence (beta selection)")
cem |
Object of class |
... |
Optional parameters. |
plot_title |
title of the graph |
Object of class CEMiTool
with beta x R squared plot
# Get example CEMiTool object data(cem) # Plot scale-free model fit as a function of the soft-thresholding beta parameter choice cem <- plot_beta_r2(cem) # Check resulting plot show_plot(cem, "beta_r2")
# Get example CEMiTool object data(cem) # Plot scale-free model fit as a function of the soft-thresholding beta parameter choice cem <- plot_beta_r2(cem) # Check resulting plot show_plot(cem, "beta_r2")
Creates a heatmap with the results of gene set enrichment analysis (GSEA) of co-expression modules
plot_gsea(cem, ...) ## S4 method for signature 'CEMiTool' plot_gsea(cem, pv_cut = 0.05)
plot_gsea(cem, ...) ## S4 method for signature 'CEMiTool' plot_gsea(cem, pv_cut = 0.05)
cem |
Object of class |
... |
Optional parameters. |
pv_cut |
P-value cut-off. Default |
Object of class CEMiTool
with GSEA plots
# Get example CEMiTool object data(cem) # Get example sample annotation file # Run GSEA on network modules cem <- mod_gsea(cem) # Plot GSEA results cem <- plot_gsea(cem) # Check resulting plot show_plot(cem, "gsea")
# Get example CEMiTool object data(cem) # Get example sample annotation file # Run GSEA on network modules cem <- mod_gsea(cem) # Plot GSEA results cem <- plot_gsea(cem) # Check resulting plot show_plot(cem, "gsea")
This function plots a histogram of the distribution of gene expression, to help assess the normality of the data.
plot_hist(cem, ...) ## S4 method for signature 'CEMiTool' plot_hist(cem, filter = FALSE)
plot_hist(cem, ...) ## S4 method for signature 'CEMiTool' plot_hist(cem, filter = FALSE)
cem |
Object of class |
... |
Optional parameters |
filter |
Logical. Whether or not to use filtered data for CEMiTool objects (Default: FALSE). |
Object of class CEMiTool
containing expression histogram
# Get example CEMiTool object data(cem) # Plot histogram cem <- plot_hist(cem) # Check results show_plot(cem, "hist")
# Get example CEMiTool object data(cem) # Plot histogram cem <- plot_hist(cem) # Check results show_plot(cem, "hist")
Creates a graph based on interactions provided
plot_interactions(cem, ...) ## S4 method for signature 'CEMiTool' plot_interactions(cem, n = 10, ...)
plot_interactions(cem, ...) ## S4 method for signature 'CEMiTool' plot_interactions(cem, n = 10, ...)
cem |
Object of class |
... |
Optional parameters. |
n |
number of nodes to label |
Object of class CEMiTool
with profile plots
# Get example CEMiTool object data(cem) # Get example gene interactions data int <- system.file("extdata", "interactions.tsv", package = "CEMiTool") int_df <- read.delim(int) # Include interaction data into CEMiTool object interactions_data(cem) <- int_df # Plot resulting networks cem <- plot_interactions(cem) # Check resulting plot show_plot(cem, "interaction")
# Get example CEMiTool object data(cem) # Get example gene interactions data int <- system.file("extdata", "interactions.tsv", package = "CEMiTool") int_df <- read.delim(int) # Include interaction data into CEMiTool object interactions_data(cem) <- int_df # Plot resulting networks cem <- plot_interactions(cem) # Check resulting plot show_plot(cem, "interaction")
Creates a graph showing the mean connectivity of genes in the network
plot_mean_k(cem, ...) ## S4 method for signature 'CEMiTool' plot_mean_k(cem, title = "Mean connectivity")
plot_mean_k(cem, ...) ## S4 method for signature 'CEMiTool' plot_mean_k(cem, title = "Mean connectivity")
cem |
Object of class |
... |
Optional parameters. |
title |
title of the graph |
Object of class CEMiTool
with connectivity plot
# Get example CEMiTool object data(cem) # Plot scale-free model fit as a function of the soft-thresholding beta parameter choice cem <- plot_mean_k(cem) # Check resulting plot show_plot(cem, "mean_k")
# Get example CEMiTool object data(cem) # Plot scale-free model fit as a function of the soft-thresholding beta parameter choice cem <- plot_mean_k(cem) # Check resulting plot show_plot(cem, "mean_k")
This plot returns a scatterplot of the mean by the variance of gene expression. A linear relationship between these values for RNAseq data suggest that an appropriate transformation such as the Variance Stabilizing Transformation should be applied.
plot_mean_var(cem, ...) ## S4 method for signature 'CEMiTool' plot_mean_var(cem, filter = FALSE)
plot_mean_var(cem, ...) ## S4 method for signature 'CEMiTool' plot_mean_var(cem, filter = FALSE)
cem |
Object of class |
... |
Optional parameters |
filter |
Logical. Whether or not to use filtered data for CEMiTool objects (Default: FALSE). |
Object of class CEMiTool
containing a mean and variance plot
# Get example CEMiTool object data(cem) # Plot mean and variance plot cem <- plot_mean_var(cem) # Check results show_plot(cem, 'mean_var')
# Get example CEMiTool object data(cem) # Plot mean and variance plot cem <- plot_mean_var(cem) # Check results show_plot(cem, 'mean_var')
Creates a bar plot with the results of module overrepresentation analysis
plot_ora(cem, ...) ## S4 method for signature 'CEMiTool' plot_ora(cem, n = 10, pv_cut = 0.05, ...)
plot_ora(cem, ...) ## S4 method for signature 'CEMiTool' plot_ora(cem, n = 10, pv_cut = 0.05, ...)
cem |
Object of class |
... |
parameters to plot_ora_single |
n |
number of modules to show |
pv_cut |
p-value significance cutoff. Default is 0.05. |
Object of class CEMiTool
with ORA plots
# Get example CEMiTool object data(cem) # Read example gmt file gmt <- read_gmt(system.file('extdata', 'pathways.gmt', package='CEMiTool')) # Run overrepresentation analysis cem <- mod_ora(cem, gmt) # Plot module gene expression profiles cem <- plot_ora(cem) # Check resulting plot show_plot(cem, "ora")
# Get example CEMiTool object data(cem) # Read example gmt file gmt <- read_gmt(system.file('extdata', 'pathways.gmt', package='CEMiTool')) # Run overrepresentation analysis cem <- mod_ora(cem, gmt) # Plot module gene expression profiles cem <- plot_ora(cem) # Check resulting plot show_plot(cem, "ora")
Creates a plot with module gene expression profiles along samples
plot_profile(cem, ...) ## S4 method for signature 'CEMiTool' plot_profile(cem, order_by_class = TRUE, center_func = "mean")
plot_profile(cem, ...) ## S4 method for signature 'CEMiTool' plot_profile(cem, order_by_class = TRUE, center_func = "mean")
cem |
Object of class |
... |
Optional parameters. |
order_by_class |
Logical. Only used if a sample
annotation file is present. Whether or not to order by the
class column in the sample annotation file (as defined by the
class_column slot in |
center_func |
Character string indicating the centrality measure to show in the plot. Either 'mean' (the default) or 'median'. |
Object of class CEMiTool
with profile plots
# Get example CEMiTool object data(cem) # Plot module gene expression profiles cem <- plot_profile(cem) # Check resulting plot show_plot(cem, "profile")
# Get example CEMiTool object data(cem) # Plot module gene expression profiles cem <- plot_profile(cem) # Check resulting plot show_plot(cem, "profile")
This function creates a normal QQ plot of the expression values.
plot_qq(cem, ...) ## S4 method for signature 'CEMiTool' plot_qq(cem, filter = FALSE)
plot_qq(cem, ...) ## S4 method for signature 'CEMiTool' plot_qq(cem, filter = FALSE)
cem |
Object of class |
... |
Optional parameters |
filter |
Logical. Whether or not to use filtered data for CEMiTool objects (Default: FALSE). |
Object of class CEMiTool
containing qqplot
# Get example CEMiTool object data(cem) # Plot quantile-quantile plot cem <- plot_qq(cem) # Check results show_plot(cem, 'qq')
# Get example CEMiTool object data(cem) # Plot quantile-quantile plot cem <- plot_qq(cem) # Check results show_plot(cem, 'qq')
Creates a dendrogram showing the similarities between samples in the expression data.
plot_sample_tree(cem, ...) ## S4 method for signature 'CEMiTool' plot_sample_tree( cem, col_vector = NULL, sample_name_column = NULL, class_column = NULL, filter = FALSE )
plot_sample_tree(cem, ...) ## S4 method for signature 'CEMiTool' plot_sample_tree( cem, col_vector = NULL, sample_name_column = NULL, class_column = NULL, filter = FALSE )
cem |
Object of class |
... |
Optional parameters. |
col_vector |
A vector of columns to use for visualizing the clustering. See Details. |
sample_name_column |
A string specifying the column to be used as sample identification. For CEMiTool objects, this will be the string specified in the sample_name_column slot. |
class_column |
A string specifying the column to be used as sample group identification. For CEMiTool objects, this will be the string specified in the class_column slot. |
filter |
Logical. Whether or not to use filtered data for CEMiTool objects (Default: FALSE). |
Object of class CEMiTool
with dendrogram or a plot object.
# Get example CEMiTool object data(cem) # Plot sample dendrogram cem <- plot_sample_tree(cem) # Check resulting plot show_plot(cem, "sample_tree")
# Get example CEMiTool object data(cem) # Plot sample dendrogram cem <- plot_sample_tree(cem) # Check resulting plot show_plot(cem, "sample_tree")
Read a GMT file
read_gmt(fname)
read_gmt(fname)
fname |
GMT file name. |
A list containing genes and description of each pathway
# Read example gmt file gmt_fname <- system.file("extdata", "pathways.gmt", package = "CEMiTool") gmt_in <- read_gmt(gmt_fname)
# Read example gmt file gmt_fname <- system.file("extdata", "pathways.gmt", package = "CEMiTool") gmt_in <- read_gmt(gmt_fname)
Modified data from a yellow fever vaccination study by Querec et al,
2009. This dataset, together with expr
can be used as input for
CEMiTool functions
data(sample_annot)
data(sample_annot)
An object of class data.frame
Querec TD, Akondy RS, Lee EK, Cao W et al. Systems biology approach predicts immunogenicity of the yellow fever vaccine in humans. Nat Immunol 2009 Jan;10(1):116-25. PMID: 19029902 PubMed
data(expr) data(sample_annot) # Run CEMiTool analysis ## Not run: cemitool(expr, sample_annot)
data(expr) data(sample_annot) # Run CEMiTool analysis ## Not run: cemitool(expr, sample_annot)
Retrive or set the sample_annotation attribute
sample_annotation(cem) ## S4 method for signature 'CEMiTool' sample_annotation(cem) sample_annotation( cem, sample_name_column = "SampleName", class_column = "Class" ) <- value ## S4 replacement method for signature 'CEMiTool' sample_annotation( cem, sample_name_column = "SampleName", class_column = "Class" ) <- value
sample_annotation(cem) ## S4 method for signature 'CEMiTool' sample_annotation(cem) sample_annotation( cem, sample_name_column = "SampleName", class_column = "Class" ) <- value ## S4 replacement method for signature 'CEMiTool' sample_annotation( cem, sample_name_column = "SampleName", class_column = "Class" ) <- value
cem |
Object of class |
sample_name_column |
A string containing the name of a column which should be used as a unique identifier for samples in the file. Only used when assigning a sample annotation data.frame. Default: "SampleName". |
class_column |
A string containing the name of a column which should be used to identify different sample groups. Only used when assigning a sample annotation data.frame. Default: "Class" |
value |
A data.frame containing the sample annotation, should have at least two columns containing the Class and the Sample Name that should match with samples in expression |
A data.frame containing characteristics of each sample.
# Get example expression data data(expr0) # Get example sample_annotation data data(sample_annot) # Initialize CEMiTool object with expression cem <- new_cem(expr0) # Add sample annotation file to CEMiTool object sample_annotation(cem, sample_name_column="SampleName", class_column="Class") <- sample_annot # Check annotation head(sample_annotation(cem))
# Get example expression data data(expr0) # Get example sample_annotation data data(sample_annot) # Initialize CEMiTool object with expression cem <- new_cem(expr0) # Add sample annotation file to CEMiTool object sample_annotation(cem, sample_name_column="SampleName", class_column="Class") <- sample_annot # Check annotation head(sample_annotation(cem))
Save plots into the directory specified by the directory
argument.
save_plots(cem, ...) ## S4 method for signature 'CEMiTool' save_plots( cem, value = c("all", "profile", "gsea", "ora", "interaction", "beta_r2", "mean_k", "sample_tree", "mean_var", "hist", "qq"), force = FALSE, directory = "./Plots" )
save_plots(cem, ...) ## S4 method for signature 'CEMiTool' save_plots( cem, value = c("all", "profile", "gsea", "ora", "interaction", "beta_r2", "mean_k", "sample_tree", "mean_var", "hist", "qq"), force = FALSE, directory = "./Plots" )
cem |
Object of class |
... |
Optional parameters One of "all", "profile", "gsea", "ora", "interaction", "beta_r2", "mean_k", "sample_tree", "mean_var", "hist", "qq". |
value |
A character string containing the name of the plot to be saved. |
force |
If the directory exists, execution will not stop. |
directory |
Directory into which the files will be saved. |
A pdf file or files with the desired plot(s)
# Get example CEMiTool object data(cem) # Plot beta x R squared graph cem <- plot_beta_r2(cem) # Save plot ## Not run: save_plots(cem, value="beta_r2", directory="./Plots")
# Get example CEMiTool object data(cem) # Plot beta x R squared graph cem <- plot_beta_r2(cem) # Save plot ## Not run: save_plots(cem, value="beta_r2", directory="./Plots")
Select genes based on variance
select_genes(expr, n_genes, filter_pval = 0.1)
select_genes(expr, n_genes, filter_pval = 0.1)
expr |
A data.frame containing expression values |
n_genes |
(Optional) Number of genes to be selected |
filter_pval |
P-value cutoff for gene selection |
A vector containing the names of selected genes
# Get example expression data data(expr0) # Filter genes expr_f <- filter_genes(expr0) # Check selected genes expr_f[1:5, 1:5] # Filter genes and apply variance stabilizing transformation expr_f2 <- filter_genes(expr0, apply_vst=TRUE) # Check results expr_f2[1:5, 1:5] # Selected genes selected <- select_genes(expr_f2) # Get data.frame with only selected genes expr_s <- expr_f2[selected, ] # Check results expr_s[1:5, 1:5]
# Get example expression data data(expr0) # Filter genes expr_f <- filter_genes(expr0) # Check selected genes expr_f[1:5, 1:5] # Filter genes and apply variance stabilizing transformation expr_f2 <- filter_genes(expr0, apply_vst=TRUE) # Check results expr_f2[1:5, 1:5] # Selected genes selected <- select_genes(expr_f2) # Get data.frame with only selected genes expr_s <- expr_f2[selected, ] # Check results expr_s[1:5, 1:5]
Retrieve CEMiTool object plots
show_plot(cem, value) ## S4 method for signature 'CEMiTool' show_plot( cem, value = c("profile", "gsea", "ora", "interaction", "beta_r2", "mean_k", "sample_tree", "mean_var", "hist", "qq") )
show_plot(cem, value) ## S4 method for signature 'CEMiTool' show_plot( cem, value = c("profile", "gsea", "ora", "interaction", "beta_r2", "mean_k", "sample_tree", "mean_var", "hist", "qq") )
cem |
Object of class |
value |
A character string containing the name of the plot to be shown. One of "profile", "gsea", "ora", "interaction", "beta_r2", "mean_k", "sample_tree", "mean_var", "hist", "qq". |
A plot corresponding to a CEMiTool analysis
# Get example CEMiTool object data(cem) # Plot beta x R squared graph cem <- plot_beta_r2(cem) # Check plot show_plot(cem, "beta_r2")
# Get example CEMiTool object data(cem) # Plot beta x R squared graph cem <- plot_beta_r2(cem) # Check plot show_plot(cem, "beta_r2")
Print a cemitool object
## S4 method for signature 'CEMiTool' show(object)
## S4 method for signature 'CEMiTool' show(object)
object |
Object of class CEMiTool |
A CEMiTool object.
Save the CEMiTool object in files
write_files(cem, ...) ## S4 method for signature 'CEMiTool' write_files(cem, directory = "./Tables", force = FALSE)
write_files(cem, ...) ## S4 method for signature 'CEMiTool' write_files(cem, directory = "./Tables", force = FALSE)
cem |
Object of class |
... |
Optional parameters |
directory |
a directory |
force |
if the directory exists the execution will not stop |
A directory containing CEMiTool results in files.
# Get example CEMiTool object data(cem) # Save CEMiTool results in files write_files(cem, directory=".", force=TRUE)
# Get example CEMiTool object data(cem) # Save CEMiTool results in files write_files(cem, directory=".", force=TRUE)