Title: | Meta-Analysis of High-Throughput Experiments Using Feature Annotations |
---|---|
Description: | Facilitates comparison of significant annotations (categories) generated on one or more feature lists. Interactive exploration is facilitated through the use of RCytoscape (heavily suggested). |
Authors: | Robert M Flight [aut, cre] |
Maintainer: | Robert M Flight <[email protected]> |
License: | MIT + file LICENSE |
Version: | 0.100.26 |
Built: | 2024-10-30 20:28:27 UTC |
Source: | https://github.com/MoseleyBioinformaticsLab/categoryCompare2 |
given the annotation_graph and a data.frame, add all of the data in the data.frame to the graph so it is available elsewhere. Note that for NA integer and numerics, the value is modified to -100, and for infinite values, it is modified to 1e100.
add_data_to_graph(graph, data)
add_data_to_graph(graph, data)
graph |
the graph to work on |
data |
the data to add to it |
graphNEL
before passing to Cytoscape, add a tooltip attribute to the graph
add_tooltip( in_graph, node_data = c("name", "description"), description, separator = "\n" )
add_tooltip( in_graph, node_data = c("name", "description"), description, separator = "\n" )
in_graph |
the graph to work with |
node_data |
which pieces of node data to use |
description |
other descriptive text to use |
separator |
what separator to use for the tooltip |
the graph with a new nodeData member "tooltip"
This class holds an annotation object that defines how annotations relate to features, as well as various pieces about each annotation
Does sensical checks when creating an annotation
object.
annotation( annotation_features, annotation_type = NULL, description = character(0), links = character(0), feature_type = NULL ) ## S4 method for signature 'annotation' show(object) annotation( annotation_features, annotation_type = NULL, description = character(0), links = character(0), feature_type = NULL )
annotation( annotation_features, annotation_type = NULL, description = character(0), links = character(0), feature_type = NULL ) ## S4 method for signature 'annotation' show(object) annotation( annotation_features, annotation_type = NULL, description = character(0), links = character(0), feature_type = NULL )
annotation_features |
list of annotation to feature relationships |
annotation_type |
a simple one word description of the annotations |
description |
character vector providing descriptive text about the annotation |
links |
character vector defining html links for each annotation (may be empty) |
feature_type |
one word description of the feature type |
object |
the annotation object |
These objects may be created by hand, or may result from specific functions to create them. Most notably, this package provides functions for creating a Gene Ontology annotation.
See the annotation
, each slot is a parameter.
annotation_features
list of annotation to feature relationships
description
character vector providing descriptive text about the annotation
counts
numeric vector of how many features are in each annotation
links
character vector defining html links for each annotation (may be empty)
annotation_type
a one word short description of the "type" of annotation
feature_type
a one word short description of the "type" of features
Given a 'categoryCompare2' annotation object, generate a JSON representation that can be used with the command line executable
annotation_2_json(annotation_obj, json_file = NULL)
annotation_2_json(annotation_obj, json_file = NULL)
annotation_obj |
the annotation object |
json_file |
the file to save it to |
the json string (invisibly)
determine the unique combinations of annotations that exist in the
significant matrix of the cc_graph
and assign each node in the graph
to a group.
determine the unique combinations of annotations that exist in the
significant matrix of the combined_statistics
and assign each
annotation to a group.
annotation_combinations(object) ## S4 method for signature 'cc_graph' annotation_combinations(object) ## S4 method for signature 'significant_annotations' annotation_combinations(object)
annotation_combinations(object) ## S4 method for signature 'cc_graph' annotation_combinations(object) ## S4 method for signature 'significant_annotations' annotation_combinations(object)
object |
the |
node_assignment
node_assignment
Creates a tabular output of annotations to genes providing lookup of which genes are contributing to a particular annotation.
annotation_gene_table( combined_enrichment, annotations = NULL, use_db = NULL, input_type = "ENTREZID", gene_info = c("SYMBOL", "GENENAME") )
annotation_gene_table( combined_enrichment, annotations = NULL, use_db = NULL, input_type = "ENTREZID", gene_info = c("SYMBOL", "GENENAME") )
combined_enrichment |
combined enrichment object |
annotations |
which annotations to grab features from |
use_db |
the annotation database |
input_type |
what type of gene id was it? |
gene_info |
what type of info to return for each gene |
data.frame
given a node_assign
, assign colors to either the independent groups
of unique annotations, or to each of the experiments independently.
assign_colors(in_assign, type = "experiment")
assign_colors(in_assign, type = "experiment")
in_assign |
the |
type |
either "group" or "experiment" |
node_assign with colors
given a cc_graph
, find communities of nodes based on their connectivity
and weights.
assign_communities(in_graph)
assign_communities(in_graph)
in_graph |
the |
list
does a binomial test
binomial_basic( positive_cases, total_cases, p_expected = 0.5, direction = "two.sided", conf_level = 0.95 )
binomial_basic( positive_cases, total_cases, p_expected = 0.5, direction = "two.sided", conf_level = 0.95 )
positive_cases |
number of positive instances |
total_cases |
total number of cases observed |
p_expected |
what is the expected probability |
direction |
which direction is the test |
conf_level |
confidence level for the confidence interval |
list
do binomial testing
binomial_feature_enrichment( binomial_features, p_expected = 0.5, direction = "two.sided", p_adjust = "BH", conf_level = 0.95, min_features = 1 )
binomial_feature_enrichment( binomial_features, p_expected = 0.5, direction = "two.sided", p_adjust = "BH", conf_level = 0.95, min_features = 1 )
binomial_features |
a binomial_features object |
p_expected |
the expected probability (default 0.5) |
direction |
which direction to do the enrichment (two.sided, less, greater) |
p_adjust |
how to correct the p-values (default is "BH") |
conf_level |
the confidence level for the confidence interval (default is 0.95) |
min_features |
a minimum number of features that are annotated to each annotation |
enriched_result
class to hold features undergoing binomail statistical testing
positivefc
the features with positive fold-changes
negativefc
the features with negative fold-changes
annotation
annotation object
the binomial results class
positivefc
the positive log-fold-changed genes, a vector of class ANY
negativefc
the negative log-fold-changed genes
annotation
list giving the annotation to feature relationship
statistics
a statistical_results
object
The categoryCompare2 package provides functions for simple enrichment and and comparison of those enrichment results.
A cc_graph
class is a graphNEL
with the added slot of
significant
, a matrix of rows (nodes / annotations) and whether
they were found to be significant in a given enrichment (columns). This
matrix is used for classifying the annotations into different groups, and
generating either pie-charts or coloring the nodes in a visualization.
constructs a cc_graph given a graphNEL
and a significant matrix.
cc_graph(graph, significant) ## S4 method for signature 'cc_graph' show(object) cc_graph(graph, significant)
cc_graph(graph, significant) ## S4 method for signature 'cc_graph' show(object) cc_graph(graph, significant)
graph |
the |
significant |
a matrix indicating which nodes are significant in which experiment |
object |
the cc_graph to show |
significant
numeric matrix of ones and zeros
For the generation of a proper annotation-annotation relationship graph, we
need to combine the annotation-feature relationships across multiple
annotation
objects
combine_annotation_features(annotation_features)
combine_annotation_features(annotation_features)
annotation_features |
list of annotation_features to combine |
list of combined annotations
Takes multiple annotation
objects and combines them so that there
is a consistent sole set for creating the cc_graph
and providing
other information about each annotation entry.
combine_annotations(annotation_list) ## S4 method for signature 'list' combine_annotations(annotation_list)
combine_annotations(annotation_list) ## S4 method for signature 'list' combine_annotations(annotation_list)
annotation_list |
one or more |
This is one of the primary workhorse functions behind categoryCompare2.
The primary function of categoryCompare
is to enable comparisons
of different enrichment analyses. To facilitate that, we must first
combine one (really, we can do this with a single) or more
enriched_result
.
combine_enrichments(...) ## S4 method for signature 'enriched_result' combine_enrichments(...) ## S4 method for signature 'list' combine_enrichments(...)
combine_enrichments(...) ## S4 method for signature 'enriched_result' combine_enrichments(...) ## S4 method for signature 'list' combine_enrichments(...)
... |
list of enriched_result |
Given lists of named character objects, and a character vector of names to be in the final object, either get the character string from the list that has the names, or check that the character string is the same across all of the lists.
combine_text(list_characters, names_out, text_id)
combine_text(list_characters, names_out, text_id)
list_characters |
list containing named character strings |
names_out |
the full list of names to use |
text_id |
what is the name for that thing being put out |
named character vector
takes an average of the overlap
and jaccard
coefficients
combined_coefficient(n1, n2)
combined_coefficient(n1, n2)
n1 |
group 1 |
n2 |
group 2 |
double
The combined_enrichment
class holds the results of combining several
enriched_result
s together, which includes the original
enriched_result
s, as well as the cc_graph
and combined annotation
objects.
enriched
list of enriched objects
annotation
annotation
where the annotation_features
have been combined across the enriched_result
statistics
combined_statistics
of both
In the case where we have a combined_enrichment
and we want
to get all of the significant annotations from each of them, and put them
together so we can start doing real meta-analysis.
combined_significant_calls(in_results, queries)
combined_significant_calls(in_results, queries)
in_results |
a |
queries |
a list of queries that can form a call object |
Note that this function returns the original combined_enrichment
object with a modified
combined_statistics
slot where the significant annotations have been added in.
combined_enrichment
object
holds the results of extracting a bunch of statistics from a combined_enrichment
into one entity. This is useful because we want to enable multiple data representations and
simple filtering on the actual data.frame
of statistics, and this provides flexibility
to enable that.
constructor function for the combined_statistics object, makes sure that empty things get initialized correctly
combined_statistics( statistic_data, which_enrichment, which_statistic, annotation_id, significant = NULL, measured = NULL, use_names = NULL ) combined_statistics( statistic_data, which_enrichment, which_statistic, annotation_id, significant = NULL, measured = NULL, use_names = NULL )
combined_statistics( statistic_data, which_enrichment, which_statistic, annotation_id, significant = NULL, measured = NULL, use_names = NULL ) combined_statistics( statistic_data, which_enrichment, which_statistic, annotation_id, significant = NULL, measured = NULL, use_names = NULL )
statistic_data |
the data.frame of statistics |
which_enrichment |
which enrichment gave the results |
which_statistic |
which statistics were calculated in each case |
annotation_id |
the annotations for which we are returning statistics |
significant |
the significant annotations |
measured |
the measured annotations |
use_names |
the order of naming |
combined_statistics
statistic_data
a data.frame
of all of the statistics from all of the enrichments
significant
a significant_annotations
object, that may be empty
which_enrichment
a vector
giving which enrichment each column of the statistics came from
which_statistic
a vector
providing which statistic each column contains
print the annotation gene table to a CSV file
csv_annotation_table(annotation_gene_table, out_file = NULL)
csv_annotation_table(annotation_gene_table, out_file = NULL)
annotation_gene_table |
list of tables |
out_file |
the file to write to |
given all the slots for an enriched_result
, checks that all
the data is self-consistent, and creates the enriched_result
object.
enriched_result(features, universe, annotation, statistics) enriched_result(features, universe, annotation, statistics)
enriched_result(features, universe, annotation, statistics) enriched_result(features, universe, annotation, statistics)
features |
the features that were differentially expressed (see details) |
universe |
all of the features that were measured |
annotation |
an |
statistics |
a |
enriched_result
features
the "features" of interest, a vector of class ANY
universe
all of the "features" in the background
annotation
list giving the annotation to feature relationship
statistics
a statistical_results
object
Takes an 'enriched_result', and converts it to the table expected by 'fgsea'. This should only be done on those that have 'gsea' as the *Enrichment Method*.
enriched_to_fgsea(in_enriched)
enriched_to_fgsea(in_enriched)
in_enriched |
the enrichment object |
data.table
Show the path to the executables, so the user can add them to whatever they want.
executable_path()
executable_path()
Extract statistical table from a single enrichment object.
extract_enrich_stats(enrichment_result)
extract_enrich_stats(enrichment_result)
enrichment_result |
the enrichment result object |
data.frame
extract all statistics for a statistical_results
object. These
can then be combined into a data.frame
that can be returned or used
to annotate the graph of annotations.
extract_statistics(in_results) ## S4 method for signature 'statistical_results' extract_statistics(in_results)
extract_statistics(in_results) ## S4 method for signature 'statistical_results' extract_statistics(in_results)
in_results |
the |
data.frame
extract all statistics from a combined_enrichment
object and
create a combined_statistics
where each statistic from the underlying
statistical_results
object in each of the enrichments
is named according to which enrichment it was in and what statistic it was.
## S4 method for signature 'combined_enrichment' extract_statistics(in_results)
## S4 method for signature 'combined_enrichment' extract_statistics(in_results)
in_results |
the |
combined_statistics
If a graph has already been generated, it may be faster to filter a previously generated one than generate a new one from significant data.
filter_annotation_graph(in_graph, comb_enrich)
filter_annotation_graph(in_graph, comb_enrich)
in_graph |
the |
comb_enrich |
the |
cc_graph
given a combined_enrichment
, generate the annotation similarity graph
generate_annotation_graph( comb_enrichment, annotation_similarity = "combined", low_cut = 5, hi_cut = 500 ) ## S4 method for signature 'combined_enrichment' generate_annotation_graph( comb_enrichment, annotation_similarity = "combined", low_cut = 5, hi_cut = 500 )
generate_annotation_graph( comb_enrichment, annotation_similarity = "combined", low_cut = 5, hi_cut = 500 ) ## S4 method for signature 'combined_enrichment' generate_annotation_graph( comb_enrichment, annotation_similarity = "combined", low_cut = 5, hi_cut = 500 )
comb_enrichment |
the combined_enrichment object |
annotation_similarity |
which similarity measure to use |
low_cut |
keep only those annotations in the graph with at least this many annotated features |
hi_cut |
keep only those annotations with less than this many annotated features |
given an annotation-feature list, generate a similarity graph between all of the annotations
generate_annotation_similarity_graph( annotation_features, similarity_type = "combined" )
generate_annotation_similarity_graph( annotation_features, similarity_type = "combined" )
annotation_features |
list where each entry is a set of features to that annotation |
similarity_type |
which type of overlap coefficient to report |
cc_graph
given a bunch of items, generate a set of colors for either single node colorings
or pie-chart annotations. Colors are generated using the hcl colorspace,
and for n_color >= 5
, the colors are re-ordered in an attempt to create
the largest contrasts between colors, as they result from being picked on a
circle in hcl space.
generate_colors(n_color)
generate_colors(n_color)
n_color |
how many colors to generate |
it often helps to have a legend displayed for reference.
generate_legend( in_assign, upper_names = TRUE, img = FALSE, width = 800, height = 400, pointsize = 70, ... )
generate_legend( in_assign, upper_names = TRUE, img = FALSE, width = 800, height = 400, pointsize = 70, ... )
in_assign |
the assign object from |
upper_names |
whether to make names uppercase for easier viewing |
img |
should a base64 encoded data uri be returned for embedding? |
width |
how wide should the image be if saving to an image |
height |
how high should it be |
pointsize |
the pointsize parameter for Cairo, determines textsize in the image |
... |
any other parameter to |
given a named vector of links, generate an actual html link formatted for output in html documents
generate_link(links)
generate_link(links)
links |
the vector of links |
character
given a group matrix and the colors for each experiment, generate the pie graphs that will be used as glyphs in Cytoscape
generate_piecharts(grp_matrix, use_color)
generate_piecharts(grp_matrix, use_color)
grp_matrix |
the group matrix |
use_color |
the colors for each experiment |
this should not be exported in the final version
list of png files that are pie graphs
given a combined_enrichment
object, get out the data.frame
either for investigation or to add data to the cc_graph
.
generate_table(comb_enrichment, link_type = "explicit") ## S4 method for signature 'combined_enrichment' generate_table(comb_enrichment, link_type = "explicit")
generate_table(comb_enrichment, link_type = "explicit") ## S4 method for signature 'combined_enrichment' generate_table(comb_enrichment, link_type = "explicit")
comb_enrichment |
the |
link_type |
should their be an "explicit" link (see details) |
the link_type
controls whether to create an "explicit" link
that is actually a column in the data.frame, or create an "implicit" html link
that is part of the @name
column in the returned data.frame. Useful
if you are embedding the data.frame in an html report.
data.frame
Generate an annotation object for genes based on an "org.*.db" object, and pulling information from it.
get_db_annotation( orgdb = "org.Hs.eg.db", features = NULL, feature_type = "ENTREZID", annotation_type = "GO" )
get_db_annotation( orgdb = "org.Hs.eg.db", features = NULL, feature_type = "ENTREZID", annotation_type = "GO" )
orgdb |
the name of the org.*.db object |
features |
which features to get annotations for |
feature_type |
which type of IDs to map (see details) |
annotation_type |
the type of annotation to grab (see details) |
This function generates a categoryCompare2
annotation object
from a Bioconductor "org.*.db" object. Even though different gene identifiers can
be used, almost all of the mappings are via ENTREZID.
The set of feature or gene keys that can be used to create the annotations include:
ENTREZID: ENTREZ gene ids
ACCNUM: genbank accession numbers
SYMBOL: gene symbols, eg ABCA1
GENENAME: gene names, eg "ATP binding cassette subfamily A member 1"
ENSEMBL: the ensembl gene ids (all start with ENSG...)
ENSEMBLPROT: ensembl protein ids (ENSP...)
ENSEMBLTRANS: ensemlb transcript ids (ENST...)
REFSEQ: reference sequence IDs, NM, NP, NR, XP, etc
UNIGENE: gene ids from UNIPROT eg Hs.88556
UNIPROT: protein ids from UNIPROT eg P80404
The set of annotations that can be mapped to features include:
GO: annotations from gene ontology
PATH: KEGG Pathway identifiers (not updated since 2011!)
CHRLOC: location on the chromosome
OMIM: mendelian inheritance in man identifiers
PMID: pubmed identifiers
PROSITE
PFAM: protein family identifiers
IPI: protein-protein interactions
For GO annotations, it is also possible to pass GO
to use all 3 sub-ontologies simultaneously,
or any combination of BP
, MF
, and CC
.
annotation object
given a statistical_results
object and some conditional expressions,
return the significant annotations
In the case where we have a combined_enrichment
and we want
to get all of the significant annotations from each of them, and put them
together so we can start doing real meta-analysis.
get_significant_annotations(in_results, ...) ## S4 method for signature 'statistical_results' get_significant_annotations(in_results, ...) ## S4 method for signature 'combined_enrichment' get_significant_annotations(in_results, ...)
get_significant_annotations(in_results, ...) ## S4 method for signature 'statistical_results' get_significant_annotations(in_results, ...) ## S4 method for signature 'combined_enrichment' get_significant_annotations(in_results, ...)
in_results |
a |
... |
conditional expressions |
Note that this function returns the original combined_enrichment
object with a modified
combined_statistics
slot where the significant annotations have been added in.
vector of significant annotation_id's
combined_enrichment
object
test_stat <- new("statistical_results", annotation_id = c("a1", "a2", "a3"), statistic_data = list(pvalues = c(a1 = 0.01, a2 = 0.5, a3 = 0.0001), counts = c(a1 = 5, a2 = 10, a3 = 1), odds = c(a1 = 20, a2 = 100, a3 = 0))) get_significant_annotations(test_stat, pvalues < 0.05) get_significant_annotations(test_stat, odds > 10) get_significant_annotations(test_stat, pvalues < 0.05, counts >= 1)
test_stat <- new("statistical_results", annotation_id = c("a1", "a2", "a3"), statistic_data = list(pvalues = c(a1 = 0.01, a2 = 0.5, a3 = 0.0001), counts = c(a1 = 5, a2 = 10, a3 = 1), odds = c(a1 = 20, a2 = 100, a3 = 0))) get_significant_annotations(test_stat, pvalues < 0.05) get_significant_annotations(test_stat, odds > 10) get_significant_annotations(test_stat, pvalues < 0.05, counts >= 1)
In the case where we have a statistical_results
and we want
to get all of the significant annotations from it
get_significant_annotations_calls(in_results, queries)
get_significant_annotations_calls(in_results, queries)
in_results |
a |
queries |
a list of queries that can form a call object |
vector of significant annotation_id's
Transforms a gocats ancestors JSON list to a GO annotation object.
gocats_to_annotation( ancestors_file = "ancestors.json", namespace_file = "namespace.json", annotation_type = "gocatsGO", feature_type = "Uniprot", feature_translation = NULL )
gocats_to_annotation( ancestors_file = "ancestors.json", namespace_file = "namespace.json", annotation_type = "gocatsGO", feature_type = "Uniprot", feature_translation = NULL )
ancestors_file |
the ancestors.json file from gocats (required) |
namespace_file |
the namespace.json file from gocats (optional) |
annotation_type |
what annotations are we making? (gocatsGO by default) |
feature_type |
what type of features are we using (assume Uniprot) |
feature_translation |
a data.frame used to convert the feature IDs |
annotation object
takes a cc_graph
object and transforms it into something that can
be visualized using visNetwork
graph_to_visnetwork( in_graph, in_assign, node_communities = NULL, use_nodes = NULL )
graph_to_visnetwork( in_graph, in_assign, node_communities = NULL, use_nodes = NULL )
in_graph |
the cc_graph object |
in_assign |
the colors generated by |
node_communities |
the communities generated by |
use_nodes |
the list of nodes to actually use |
list
Performs gene-set enrichment analysis using the 'fgsea' package.
gsea_feature_enrichment( gsea_features, min_features = 15, max_features = 500, return_type = "cc2", ... )
gsea_feature_enrichment( gsea_features, min_features = 15, max_features = 500, return_type = "cc2", ... )
gsea_features |
a GSEA features object |
min_features |
the minimum number of features for an annotation (default = 15) |
max_features |
the maximum number of features for an annotation (default = 500) |
return_type |
what type of object should be returned? ("cc2" or "fgsea") |
... |
other 'fgsea' options |
The runtime is dependent on the maximum size of the provided annotation, so the authors of 'fgsea' recommend a maximum size of 500. In addition, to calculate statistics, a minimum size of annotated features are required. Going below 15 may not be advised. If you want to use other 'fgsea' functions, it is recommended to set 'return_type = "fgsea"'. Otherwise, you should keep the default of "cc2".
enriched_result
fgsea::fgsea
class to hold features undergoing GSEA
ranks
a named vector of ranks
annotation
annotation object
class to hold features undergoing hypergeometric enrichment
significant
the significant features
universe
all of the features measured
annotation
annotation object
does a hypergeometric enrichment test
hypergeometric_basic( num_white, num_black, num_drawn, num_white_drawn, direction = "over" )
hypergeometric_basic( num_white, num_black, num_drawn, num_white_drawn, direction = "over" )
num_white |
number of white balls in urn |
num_black |
number of black balls in urn |
num_drawn |
number of balls taken from urn |
num_white_drawn |
number of white balls taken from urn |
direction |
which direction is the test |
list
do hypergeometric enrichment
hypergeometric_feature_enrichment( hypergeom_features, direction = "over", p_adjust = "BH", min_features = 1 )
hypergeometric_feature_enrichment( hypergeom_features, direction = "over", p_adjust = "BH", min_features = 1 )
hypergeom_features |
a hypergeometric_features object |
direction |
which direction to do the enrichment (over or under) |
p_adjust |
how to correct the p-values (default is "BH") |
min_features |
how many features should be annotated before testing it? |
The min_features
argument here applies to the minumum number of features an annotation has
from the universe of features supplied, not the minumum number of features from the differential
list.
For more about the p-value adjustment, see stats::p.adjust
enriched_result
move executables to user location, default is ~/bin and changes their permissions to make them executable.
install_executables(path = "~/bin")
install_executables(path = "~/bin")
path |
the path to put the executable scripts |
the listing of the files.
calculates similarity of two groups of objects using "jaccard" coefficient, defined as:
jaccard_coefficient(n1, n2)
jaccard_coefficient(n1, n2)
n1 |
group 1 |
n2 |
group 2 |
length(intersect(n1, n2)) / min(c(length(n1), length(n2)))
double
Given a JSON based annotation object, read it in and create the 'annotation' for actually doing enrichment.
json_2_annotation(json_file)
json_2_annotation(json_file)
json_file |
the json annotation file |
annotation object
Given a JSON file of features to annotations, reverse to turn it into annotations to features, and optionally add some meta-information about them.
json_annotation_reversal( json_file, out_file = "annotations.json", feature_type = NULL, annotation_type = NULL )
json_annotation_reversal( json_file, out_file = "annotations.json", feature_type = NULL, annotation_type = NULL )
json_file |
the json file to use |
out_file |
the json file to write out to |
feature_type |
the type of features |
annotation_type |
the type of annotations |
the json object, invisibly
print the annotation gene table in knitr::kable
format
kable_annotation_table(annotation_gene_table, header_level = 3, cat = TRUE)
kable_annotation_table(annotation_gene_table, header_level = 3, cat = TRUE)
annotation_gene_table |
list of tables |
header_level |
what header level should the labels be done at? |
cat |
whether to write it directly, or just return the table for later |
character
Determine the label of a community based on the most generic member of each community, which is defined as being the one with the most annotations.
label_communities(community_defs, annotation)
label_communities(community_defs, annotation)
community_defs |
the communities from |
annotation |
the annotation object used for enrichment |
list
Provided a list, and a condition, returns the logical indices into the named
part of the list provided. Uses subset
like non-standard evaluation
so that we can define appropriate expressions.
multi_query_list(list_to_query, ...)
multi_query_list(list_to_query, ...)
list_to_query |
the list to run the query on |
... |
the expressions that do the queries |
logical "&" of all queries
The node_assign
class holds the unique annotation combinations and the
assignment of the nodes to those combinations for use in visualization.
groups
the unique groups, as a logical matrix
assignments
named character vector providing association with groups
description
named character vector providing a description to group
colors
named character vector of hex colors for groups or experiments
color_type
whether doing group or experiment based colors
pie_locs
if doing experiment colors, then pie graphs were generated here
calculates the similarity using the "overlap" coefficient, which is
overlap_coefficient(n1, n2)
overlap_coefficient(n1, n2)
n1 |
group 1 of objects |
n2 |
group 2 of objects |
length(intersect(n1, n2)) / length(union(n1, n2))
double
given a RCy3 network connection, remove edges according to provided values.
remove_edges(edge_obj, cutoff, edge_attr = "weight", value_direction = "under") ## S4 method for signature 'character,numeric' remove_edges(edge_obj, cutoff, edge_attr = "weight", value_direction = "under") ## S4 method for signature 'cc_graph,numeric' remove_edges(edge_obj, cutoff, edge_attr = "weight", value_direction = "under")
remove_edges(edge_obj, cutoff, edge_attr = "weight", value_direction = "under") ## S4 method for signature 'character,numeric' remove_edges(edge_obj, cutoff, edge_attr = "weight", value_direction = "under") ## S4 method for signature 'cc_graph,numeric' remove_edges(edge_obj, cutoff, edge_attr = "weight", value_direction = "under")
edge_obj |
cc_graph |
cutoff |
the cutoff to use |
edge_attr |
which attribute to use |
value_direction |
remove edges with value under or over |
nothing
cc_graph
show binomial_result
## S4 method for signature 'binomial_result' show(object)
## S4 method for signature 'binomial_result' show(object)
object |
the binomial_result object to show |
show combined_statistics
## S4 method for signature 'combined_statistics' show(object)
## S4 method for signature 'combined_statistics' show(object)
object |
show enriched_result
## S4 method for signature 'enriched_result' show(object)
## S4 method for signature 'enriched_result' show(object)
object |
the enriched_result object to show |
show node_assign
## S4 method for signature 'node_assign' show(object)
## S4 method for signature 'node_assign' show(object)
object |
the node_assign to see |
show signficant_annotations
## S4 method for signature 'significant_annotations' show(object)
## S4 method for signature 'significant_annotations' show(object)
object |
the significant annotations object to show |
The significant_annotations
class holds which annotations from which
enrichment were both measured and significant. Each of these
slots is a logical matrix with rows named by annotation_id and
columns named by the names of the enriched_result
that was combined.
Makes a new significant_annotation while checking that everything is valid.
significant_annotations(significant, measured, sig_calls = NULL) significant_annotations(significant, measured, sig_calls = NULL)
significant_annotations(significant, measured, sig_calls = NULL) significant_annotations(significant, measured, sig_calls = NULL)
significant |
logical matrix of annotations (rows) and experiments (columns) |
measured |
logical matrix of annotations (rows) and experiments (columns) |
sig_calls |
character vector of deparsed calls that resulted in signficant and measured |
significant
logical matrix
measured
logical matrix
sig_calls
character representations of calls used to filter the data
This class holds the part of an enrichment that is the statistical results.
It has two pieces, a list
of statistics
that is a named list
with the actual numerical results of applying the statistics. The other piece
is the annotation_id
vector defining which entry in each vector of the
statistics
is.
statistic_data
list of numerical statistics
annotation_id
vector of ids
method
how the statistics were calculated
Creates a table from the annotation graph, and if provided, adds the community information to the table.
table_from_graph(in_graph, in_assign = NULL, community_info = NULL)
table_from_graph(in_graph, in_assign = NULL, community_info = NULL)
in_graph |
the |
in_assign |
the |
community_info |
the |
data.frame
given a graph, and the node assignments, visualize the graph in cytoscape for manipulation
vis_in_cytoscape(in_graph, in_assign, description = "cc2 enrichment")
vis_in_cytoscape(in_graph, in_assign, description = "cc2 enrichment")
in_graph |
the cc_graph to visualize |
in_assign |
the node_assign generated |
description |
something descriptive about the vis (useful when lots of different visualizations) |
something
Visualize a cc_graph
in visNetwork
, with selection for communities
if that exists.
vis_visnetwork(in_graph_info)
vis_visnetwork(in_graph_info)
in_graph_info |
the graph structure from |