vignettes/v5_term_similarity.Rmd
v5_term_similarity.Rmd
The methods of semantic similarity implemented in simona are mainly from the supplementary file of the paper “Mazandu et al., Gene Ontology semantic similarity tools: survey on features and challenges for biological knowledge discovery. Briefings in Bioinformatics 2017”. Original denotations have been slightly modified to make them more consistent. Also more explanations have been added in this vignette.
The following denotations will be used throughout the vignette. The denotations are mainly from Mazandu 2017 only with small modifications.
Denotation | Description |
---|---|
r | The root term of the DAG. In simona there is always one root term. |
δ(x) | The depth of a term x in the DAG, which is the longest distance from root r. |
δs(x) | The length of the longest path from root r to a term x via term s. |
δmax | The maximal depth in the DAG. |
η(x) | The height of term x in the DAG, which is the longest finite distance to leaf terms. |
Cs | The set of child terms of term s. |
Ps | The set of parent terms of term s. |
As | The set of ancestor terms of term s. |
A+s | The set of ancestor terms of term s, including s itself. |
Ds | The set of offspring terms of term s. |
D+s | The set of offspring terms of term s, including s itself. |
Ls | The set of leaf terms that term s can reach. |
|A| | Number of elements in set A. |
Dsp(a,b) | The shortest distance bewteen a and b. |
len(a,b) | The longest distance bewteen a and b. |
lens(a,b) | The length of the longest path from a and b via s. |
CA(a,b) | The set of common ancestors of term a and b, i.e. CA(a,b)=A+a∩A+b. |
LCA(a,b) | Lowest common ancestor of a and b, which is the common ancestor with the largest depth in DAG, i.e. argmaxt∈CA(a,b)δ(t) There might be more than one LCA terms for given two terms, to simplify the calculation, the one with the longest distance (the default) to a and b is used. |
NCA(a,b) | Nearest common ancestor of a and b, i.e. argmint∈CA(a,b)(Dsp(t,a)+Dsp(t,b)) If there are more than one NCA terms, the one with the largest depth (the lowest one) is used. |
MICA(a,b) | Most informative common ancestor of a and b, i.e. argmaxt∈CA(a,b)(IC(t)) There might be more than one MICA terms for given two terms, the one with the longest distance (the default) to a and b is used. |
Gs | The set of annotated items on term s. |
Assume term a is an ancestor of term b, Dsp(a,b) (the order of a and b does not matter) is the normal shortest distance from a to b in a directed graph. The definition is similar for len(a,b).
If term a and b are not in offspring/ancestor relationship, i.e. a is not an ancestor of b, and b is not an ancestor of a, then
Dsp(a,b)=mint∈CA(a,b)(Dsp(t,a)+Dsp(t,b))len(a,b)=maxt∈CA(a,b)(len(t,a)+len(t,b))
The wrapper function term_sim()
calculates semantic similarities between terms in the DAG with a specific method. Note the method name can be partially matched. control
argument controls parameters for specific methods.
term_sim(dag, terms, method = ..., control = list(...))
All supported term similarity methods are:
## [1] "Sim_Lin_1998" "Sim_Resnik_1999" "Sim_FaITH_2010"
## [4] "Sim_Relevance_2006" "Sim_SimIC_2010" "Sim_XGraSM_2013"
## [7] "Sim_EISI_2015" "Sim_AIC_2014" "Sim_Zhang_2006"
## [10] "Sim_universal" "Sim_Wang_2007" "Sim_GOGO_2018"
## [13] "Sim_Rada_1989" "Sim_Resnik_edge_2005" "Sim_Leocock_1998"
## [16] "Sim_WP_1994" "Sim_Slimani_2006" "Sim_Shenoy_2012"
## [19] "Sim_Pekar_2002" "Sim_Stojanovic_2001" "Sim_Wang_edge_2012"
## [22] "Sim_Zhong_2002" "Sim_AlMubaid_2006" "Sim_Li_2003"
## [25] "Sim_RSS_2013" "Sim_HRSS_2013" "Sim_Shen_2010"
## [28] "Sim_SSDD_2013" "Sim_Jiang_1997" "Sim_Kappa"
## [31] "Sim_Jaccard" "Sim_Dice" "Sim_Overlap"
## [34] "Sim_Ancestor"
This type of methods consider a special ancestor term c of terms a and b, which has the highest IC among all a and b’s ancestor terms. Term c is called the most informative common ancestor (MICA) which can be given by:
IC(c)=maxt∈A+a∩A+bIC(t)
So if two terms are identical, MICA is the term itself, and if two terms have ancestor/offspring relationship, MICA is the ancestor term.
In the following sections, if not specially mentioned, c is always referred to the MICA of a and b.
The similarity is calculated as the IC of the MICA term c normalized by the average of the IC of the two terms:
Sim(a,b)=IC(c)(IC(a)+IC(b))/2=2∗IC(c)IC(a)+IC(b)
term_sim(dag, terms, method = "Sim_Lin_1998")
Paper link: https://dl.acm.org/doi/10.5555/645527.657297.
IC of the MICA term itself IC(c) can be a measure of how similar two terms are, but its range is not in [0, 1]
. There are several ways to normalize IC(c) to the range of [0, 1]
. Note some of the normalization methods are restricted to IC_annotation as the IC method.
It is normalized to the possible maximal IC value where a term only has one item annotated.
Sim(a,b)=IC(c)−log(1/N)=IC(c)logN
where N is the total number of items annotated to the whole DAG.
It is similar to Nunif, but normalized to the maximal IC of all terms in the DAG. If there is a term with only one item annotated, Nmax is identical to the Nunif method.
Sim(a,b)=IC(c)ICmax
IC(c) is normalized by the maximal IC of term a and b.
Sim(a,b)=IC(c)max{IC(a),IC(b)}
Paper link: https://doi.org/10.1613/jair.514, https://doi.org/10.1186/1471-2105-9-S5-S4, https://doi.org/10.1186/1471-2105-11-562, https://doi.org/10.1155/2013/292063.
The normalization method can be set with the norm_method
parameter:
term_sim(dag, terms, method = "Sim_Resnik_1999", control = list(norm_method = "Nunif")) term_sim(dag, terms, method = "Sim_Resnik_1999", control = list(norm_method = "Nmax")) term_sim(dag, terms, method = "Sim_Resnik_1999", control = list(norm_method = "Nunivers"))
It is calculated as:
Sim(a,b)=IC(c)IC(a)+IC(b)−IC(c)
The relation between the FaITH_2010 similarity and Lin_1998 similarity is:
SimFaITH=SimLin2−SimLin
term_sim(dag, terms, method = "Sim_FaITH_2010")
Paper link: https://doi.org/10.1007/978-3-642-17746-0_39.
If thinking Lin_1998 is a measure of how close term a and b are to their MICA term c, the relevance method corrects it by multiplying a factor which considers the specificity of how c brings the information. The factor is calculated as 1−p(c) where p(c) is the annotation-based probability p(c)=k/N where k is the number of items annotated to c and N is the total number of items annotated to the DAG. Then under the Relevance method, the corrected IC of c is:
ICcorrected(c)=(1−p(c))∗IC(c)
If using Lin_1998 as the similarity method, the corrected version Relevance similarity is:
Sim(a,b)=2∗ICcorrected(c)IC(a)+IC(b)=(1−p(c))∗2∗IC(c)IC(a)+IC(b)=(1−p(c))∗SimLin(a,b)
The term p(c) requires that terms should be annotated to items. However, it can be extended to more general scenarios:
ICcorrected(c)=(1−exp(−IC(x)))∗IC(c)
term_sim(dag, terms, method = "Sim_Relevance_2006")
Paper link: https://doi.org/10.1186/1471-2105-7-302
The SimIC_2010 method is an improved correction method of the Relevance method because the latter works badly when p(c) is very small. E.g., when 1−p(c) is used as a correction factor, it cannot nicely distinguish p(c)=0.01 and p(c)=0.001 because for both 1−p(c) are very close to 1.
The SimIC_2010 correction factor for MICA term c is:
ICcorrected(c)=1−11−log(p(c))∗IC(c)
Then the similarity (if using Lin_1998 as the original similarity method) is:
Sim(a,b)=(1−11−log(p(c)))∗SimLin(a,b)
Similarly, it can be generalized to:
Sim(a,b)=IC(x)1+IC(x)∗SimLin(a,b)
term_sim(dag, terms, method = "Sim_SimIC_2010")
Paper link: https://doi.org/10.48550/arXiv.1001.0958.
Being different from the Relevance and SimIC_2010 methods that only use the IC of the MICA term, the XGraSM_2013 method as well as the next method use IC of a subset of common ancestor terms of a and b, and it uses the mean IC of them. The subset of common ancestor may have different names for different methods.
XGraSM_2013 is the simplest one which uses informative common ancestors (ICA) where IC of the common ancestor is positive.
ICA(a,b)={c∈A+a∩A+b:IC(c)>0}
And mean IC among all ICA terms:
ICmean=1|ICA(a,b)|∑t∈ICA(a,b)IC(t)
And applying Lin_1998 method, the semantic similarity is:
Sim(a,b)=2∗ICmeanIC(a)+IC(b)
term_sim(dag, terms, method = "Sim_XGraSM_2013")
Paper link: https://doi.org/10.1186/1471-2105-14-284
It selects a specific subset of common ancestors of terms a and b. It only selects a common ancestor c which can reach a or b via one of its child terms that does not belong to the common ancestors (mutual exclusively in a’s ancestors or in b’s ancestors). The set of the selected common ancestors is called the exclusively inherited common ancestors (EICA).
EICA(a,b)={c∈Aa∩Ab:Cc∩((Aa∪Ab)−(Aa∩Ab)≠∅)}
And mean IC among all EICA terms:
ICmean=1|EICA(a,b)|∑t∈EICA(a,b)IC(t)
And applying Lin_1998 method, the semantic similarit is:
Sim(a,b)=2∗ICmeanIC(a)+IC(b)
term_sim(dag, terms, method = "Sim_EISI_2015")
Paper link: https://doi.org/10.1016/j.gene.2014.12.062
It uses the aggregate information content from ancestors. First define the semantic weight denoted as Sw of a term t in the DAG:
Sw(t)=11+exp(−1IC(t))
Then the similarity is calculated as the fraction of aggegation from common ancestors and the average aggregation from ancestors of a and b separately.
Sim(a,b)=2∗∑t∈A+a∩A+bSw(t)∑t∈A+aSw(t)+∑t∈A+bSw(t)
term_sim(dag, terms, method = "Sim_AIC_2014")
Paper link: https://doi.org/10.1109/tcbb.2013.176.
It uses the IC_Zhang_2006 IC method and uses Lin_1998 similarity method to calculate similarities:
Sim(a,b)=2∗ICZhang(c)ICZhang(a)+ICZhang(b)
term_sim(dag, terms, method = "Sim_Zhang_2006")
It uses the IC_universal IC method and uses the Nunivers method to calculate similarities:
Sim(a,b)=2∗ICUnivers(c)max{ICUnivers(a),ICUnivers(b)}
term_sim(dag, terms, method = "Sim_universal")
Similar as the Sim_AIC_2014 method, it is also aggregation from ancestors, but it uses the “S-value” introduced in the IC_Wang_2007 sectionn in 4. Information content.
Sim(a,b)=∑t∈A+a∩A+b(Sa(t)+Sb(t))∑t∈A+aSa(t)+∑t∈A+bSb(t)
The contribution of different semantic relations can be set with the contribution_factor
parameter. The value should be a named numeric vector where names should cover the relations defined in relations set in create_ontology_DAG()
. For example, if there are two relations “relation_a” and “relation_b” set in the DAG, the value for contribution_factor can be set as:
term_sim(dag, terms, method = "Sim_Wang_2007", control = list(contribution_factor = c("relation_a" = 0.8, "relation_b" = 0.6)))
By default 0.8 is set for “is_a” and 0.6 for “part_of”.
If you are not sure what types of relations have been set, simply type the dag
object. The relation types will be printed there.
Paper link: https://doi.org/10.1093/bioinformatics/btm087.
It is very similar as Sim_Wang_2007 except there is a correction for the contribution factor. When calculating the “S-value” introduced in the IC_Wang_2007 sectionn in 4. Information content, for a parent and a child, the weight variable we is directly determined by the relation type, e.g, 0.8 for “is_a”. In Sim_GOGO_2018, the number of child terms is also considered for we:
we=1c+|Ct|+w0
where |Ct| is the number of child terms of the parent t, w0 is the original contribution factor directly assigned for each relation type. c is selected to ensure we≤1 (assuming minimal number of children is 1), which is normally:
c=max{w0}1−max{w0}
By default, 0.4 is assigned for “is_a” and 0.3 is assigned for “part_of”, c is set to 2/3 (solve 1 = 1/(c + 1) + 0.4
).
term_sim(dag, terms, method = "Sim_GOGO_2018", control = list(contribution_factor = c("relation_a" = 0.4, "relation_b" = 0.3)))
Paper link: https://doi.org/10.1038/s41598-018-33219-y.
Methods introduced in this section relies on the distance between terms. Many methods are defined originally based on the shortest distance between two terms. This section extends them to also support their longest distance via the LCA term.
It is based on the distance between term a and b. It is defined as:
Sim(a,b)=11+Dsp(a,b)
which is based on the shortest distance between a and b. Optionally, the distance can also be the longest distance via the LCA term c.
Sim(a,b)=11+lenc(a,b)=11+len(c,a)+len(c,b)
There is a parameter distance which takes value of "longest_distances_via_LCA"
(the default) or "shortest_distances_via_NCA"
:
term_sim(dag, terms, method = "Sim_Rada_1989", control = list(distance = "shortest_distances_via_NCA"))
Paper link: https://doi.org/10.1109/21.24528.
It is a normalized distance:
Sim(a,b)=1−Dsp(a,b)2∗δmax
where 2∗δmax can be thought as the possible maximal distance between two terms in the DAG.
Similarly, the distance can also be the longest distance via LCA, then it is consistent with the definition of δmax which are both based on the longest distance.
Sim(a,b)=1−lenc(a,b)2∗δmax
There is a parameter distance which takes value of "longest_distances_via_LCA"
(the default) or "shortest_distances_via_NCA"
:
term_sim(dag, terms, method = "Sim_Resnik_edge_2005", control = list(distance = "shortest_distances_via_NCA"))
Paper link: https://doi.org/10.1145/1097047.1097051.
It is similar as the Sim_Resnik_edge_2005 method, but it applies log-transformation on the distance and the depth:
Sim(a,b)=1−log(Dsp(a,b))log(2∗δmax)
where 2∗δmax can be thought as the possible maximal distance between two terms in the DAG.
Similarly, the distance can also be the longest distance via LCA, then it is consistent with the definition of δmax which are both based on the longest distance.
Sim(a,b)=1−log(lenc(a,b))log(2∗δmax)
There is a parameter distance which takes value of "longest_distances_via_LCA"
(the default) or "shortest_distances_via_NCA"
:
term_sim(dag, terms, method = "Sim_Leocock_1998", control = list(distance = "shortest_distances_via_NCA"))
Paper link: https://ieeexplore.ieee.org/document/6287675.
It is based on the depth of the LCA term c and the longest distance between term a and b via c:
Sim(a,b)=2∗δ(c)len(c,a)+len(c,b)+2∗δ(c)=2∗δ(c)lenc(a,b)+2∗δ(c)
And it can also be written in the Lin_1998 form:
Sim(a,b)=2∗δ(c)δ(c)+len(c,a)+δ(c)+len(c,b)=2∗δ(c)δc(a)+δc(b)
where in the denominator are the depths of a and b via c.
term_sim(dag, terms, method = "Sim_WP_1994")
Paper link: https://doi.org/10.3115/981732.981751.
It is a correction of the Sim_WP_1994 method. The correction factor for term a and b regarding to their LCA term c is:
Sim(a,b)=CF(a,b)∗SimWP(a,b)
The correction factor CF(a,b) is based whether a and b are in ancestor/offspring relationship or not.
CF(a,b)={min{δ(a),δ(b)}−δ(c)=min{len(c,a),len(c,b)}a and b are not ancestor-offspring11+|δ(a)−δ(b)|=11+len(a,b)a and b are ancestor-offspring
term_sim(dag, terms, method = "Sim_Slimani_2006")
Paper link: https://zenodo.org/record/1075130.
It is also a correction of the Sim_WP_1994 method. The correction factor for term a and b is:
CF(a,b)={1a and b are not ancestor-offspringexp(−Dsp(a,b)δmax))a and b are ancestor-offspring
Dsp can be replaced with len(a,b) if the longest distance is used.
There is a parameter distance which takes value of "longest_distances_via_LCA"
(the default) or "shortest_distances_via_NCA"
:
term_sim(dag, terms, method = "Sim_Shenoy_2012", control = list(distance = "shortest_distances_via_NCA"))
Paper link: https://doi.org/10.48550/arXiv.1211.4709.
It is very similar to the Sim_WP_1994 method:
Sim(a,b)=δ(c)len(c,a)+len(c,b)+δ(c)=δ(c)δ(c)+len(c,a)+δ(c)+len(c,b)−δ(c)=δ(c)δc(a)+δc(b)−δ(c)
And the relationship to SimWP is:
SimPekar(a,b)=SimWP(a,b)2−SimWP(a,b)
term_sim(dag, terms, method = "Sim_Pekar_2002")
Paper link: https://aclanthology.org/C02-1090/.
It is purely based on the depth of term a, b and their LCA term c.
Sim(a,b)=δ(c)δ(a)+δ(b)−δ(c)
The similarity value might be negative because there is no restrction that the path from root to a or b must pass c.
term_sim(dag, terms, method = "Sim_Stojanovic_2001")
Paper link: https://doi.org/10.1145/500737.500762.
It is calculated as:
Sim(a,b)=len(r,c)2lenc(r,a)∗lenc(r,b)=δ(c)2δc(a)∗δc(b)
where r is the root term.
term_sim(dag, terms, method = "Sim_Wang_edge_2012")
Paper link: https://doi.org/10.1186/1477-5956-10-s1-s18.
For a term x, it first calculates a “mile-stone” value based on the depth as
m(x)=2−δ(x)−1
The the distance bewteen term a and b via LCA term c is:
D(a,b)=D(c,a)+D(c,b)=m(c)−m(a)+m(c)+m(b)=2−δ(c)−2−δ(a)−1−2−δ(b)−1
We can change original δ(a) and δ(b) to δc(a) and δc(b) to require that the depth to reach a and b should go through c. Then above equation becomes
D(a,b)=2−δ(c)−2−δc(a)−1−2−δc(b)−1=2−δ(c)−2−δ(c)−len(c,a)−1−2−δ(c)−len(c,b)−1=2−δ(c)(1−2−len(c,a)−1−2−len(c,b)−1)
Then when a=b (the two terms are identical), D(a,b)=0 and when c=r (common ancestor only includes root) and len(r,a)→∞, len(r,b)→∞ (root has infinite distance to the terms), D(a,b) reaches maximal of 1. So the similarity
Sim(a,b)=1−D(a,b)
ranges between 0 and 1.
term_sim(dag, terms, method = "Sim_Zhong_2002")
Paper link: https://doi.org/10.1007/3-540-45483-7_8.
It also takes accout of the distance between term a and b, as well as the depth of the LCA term c in the DAG. The distance is calculated as:
D(a,b)=log(1+Dsp(a,b)∗(σmax−σ(c)))
To scale D(a,b) into the range of [0, 1]
, we can calculate the smallest value as zero when a=b. D(a,b) reaches maximal when Dsp(a,b) reach possible maximal which is 2∗δmax. Then we can define the maximal of D(a,b) as
Dmax=log(1+2∗δmax∗δmax)
And the similarity is:
Sim(a,b)=1−D(a,b)/Dmax
There is a parameter distance which takes value of "longest_distances_via_LCA"
(the default) or "shortest_distances_via_NCA"
:
term_sim(dag, terms, method = "Sim_AlMubaid_2006", control = list(distance = "shortest_distances_via_NCA"))
Paper link: https://doi.org/10.1109/IEMBS.2006.259235.
It is similar to the Sim_AlMubaid_2006 method, but uses a non-linear form:
Sim(a,b)=exp(−0.2∗Dsp(a,b))∗tanh(0.6∗δ(c))
There is a parameter distance which takes value of "longest_distances_via_LCA"
(the default) or "shortest_distances_via_NCA"
:
term_sim(dag, terms, method = "Sim_Li_2003", control = list(distance = "shortest_distances_via_NCA"))
Paper link: https://doi.org/10.1109/TKDE.2003.1209005.
Hybrid methods use both DAG structure information and IC.
The similarity is adjusted by the positions of term a, b and the LCA term c in the DAG. The similarity is defined as:
Sim(a,b)=δmaxδmax+Dsp(a,b)∗αα+β
where Dsp(a,b) can also be the longest distance via LCA. α and β in the second term are defined as:
α=δ(c)β=min{η(a),η(b)}
where α is the depth of LCA, β corresponds to the distance to leaves, which is the smaller height of a and b in the DAG.
There is a parameter distance which takes value of "longest_distances_via_LCA"
(the default) or "shortest_distances_via_NCA"
:
term_sim(dag, terms, method = "Sim_RSS_2013", control = list(distance = "shortest_distances_via_NCA"))
Paper link: https://doi.org/10.1371/journal.pone.0066745.
It is similar to the Sim_RSS_2013 method, but it uses information content instead of the distance to adjust the similarity.
It first defines the semantic distance between term a and b as the sum of the distance to their MICA term c:
D(a,b)=D(c,a)+D(c,b)
And the distance between an ancestor to a term is:
D(c,a)=IC(a)−IC(c)D(a,b)=D(c,a)+D(c,b)=IC(a)+IC(b)−2∗IC(c)
Similarly, the similarity is also corrected by the position of MICA term and a, b in the DAG:
Sim(a,b)=11+D(a,b)∗αα+β
where
α=IC(c)
And beta is the average of the maximal semantic distance of a and b to leaves.
β=D(a,la)+D(b,lb)2=IC(la)−IC(a)+IC(lb)−IC(b)2
where la or lb is the leaf with the highest IC that a or b can reach (i.e. the most informative leaf)
IC(la)=maxz∈L(a)IC(z)
term_sim(dag, terms, method = "Sim_HRSS_2013")
Paper link: https://doi.org/10.1371/journal.pone.0066745.
It is based on the information content of terms on the path connecting term a and b via their MICA term c.
Denote a list of terms a, ..., c, ..., b
which are composed by the shortest path from c to a and from c to b, the distance between a and b is the sum of 1/IC of the terms on the path. Denote Lc(a,b) as the set of terms on the shortest path connecting a and b via the MICA term c, the similarity is:
Sim(a,b)=1−arctan(∑x∈Lc(a,b)1IC(x))π/2
The path Lc(a,b) can also be defined as the longest path via MICA. The distance
parameter controls which type of paths to use.
term_sim(dag, terms, method = "Sim_Shen_2010", control = list(distance = "longest_distances_via_LCA"))
Paper link: https://doi.org/10.1109/BIBM.2010.5706623.
It is similar to the Sim_Shen_2010 method which also sums information along the path passing through the LCA term. Instead of summing the information contents, the Sim_SSDD_2013 method sums up a so-called “T-value” which relies on the DAG structure.
Denote Lc(a,b) as the set of terms on the shortest path connecting a and b via the LCA term c, the similarity is calculated as:
Sim(a,b)=1−arctan(∑x∈Lc(a,b)T(x))π/2
The T-value T(x) depends on the DAG structure which considers both parents and children of x. The definition of T(x) is:
T(x)={1if x is a root1|Px|∑t∈Px(w∗T(t))otherwise
which means T-value of a term is an average of the weighted T-values of its parents. The weight w measures the fraction of information a parent t transmitting to downstream of the DAG via x, defined as:
w=|D+x||D+t|
w≤1 as all offsprings of x are also offspring of its parent t.
The path Lc(a,b) can also be defined as the longest path via MICA. The distance
parameter controls which type of paths to use.
term_sim(dag, terms, method = "Sim_SSDD_2013", control = list(distance = "longest_distances_via_LCA"))
Paper link: https://doi.org/10.1016/j.ygeno.2013.04.010.
First semantic distance between term a and b via MICA term c is defined as:
D(a,b)=IC(a)+IC(b)−2∗IC(c)
Then there are several normalization methods to change the distance to similarity and to scale it into the range of [0, 1]
.
"max"
: 1−D(a,b)2∗ICmax
"Couto"
: min{1,D(a,b)ICmax}
"Lin"
: 1−D(a,b)IC(a)+IC(b) which is the same as the Sim_Lin_1998 method"Garla"
: 1−log(D(a,b)+1)log(2∗ICmax+1)
"log-Lin"
: 1−log(D(a,b)+1)log(IC(a)+IC(b)+1)
"Rada"
: 11+D(a,b)
The normalization methods can be set via the parameter norm_method
:
term_sim(dag, terms, method = "Sim_Jiang_1997", control = list(norm_method = "max")) term_sim(dag, terms, method = "Sim_Jiang_1997", control = list(norm_method = "Couto")) term_sim(dag, terms, method = "Sim_Jiang_1997", control = list(norm_method = "Lin")) term_sim(dag, terms, method = "Sim_Jiang_1997", control = list(norm_method = "Garla")) term_sim(dag, terms, method = "Sim_Jiang_1997", control = list(norm_method = "log-Lin")) term_sim(dag, terms, method = "Sim_Jiang_1997", control = list(norm_method = "Rada"))
Paper link: https://aclanthology.org/O97-1002/.
Denote A and B as the sets of items annotated to term a and b, and U as the universe set of all items annotated to the DAG.
The definition of kappa coeffient is a little bit complex. First let’s format the two sets into a contigency table:
In set B | |||
Yes | No | ||
In set A | Yes | a | b |
No | c | d |
where a, b, c, d are the numbers of items that fall in each category.
Let’s calculate pobs (probability of observed agreement, both yes or both no) and pexp (probability of expected agreement) as:
pobs=a+da+b+c+dpYes=a+ba+b+c+d∗a+ca+b+c+dpNo=c+da+b+c+d∗b+da+b+c+dpexp=pYes+pNo
where pobs is the probability of an item in both sets or neither in both sets, pYes is the probability of an item in both sets by random (by assuming the events of an item in set A and set B are independent), pNo is the probability of an item not in the two sets by random, and pexp is the probability of an item either both in the two sets or not in the two sets by random.
The kappa coeffcient is calculated as:
Sim(a,b)=Kappa(a,b)=pobs−pexp1−pexp
Note the Kappa coeffcient is possible to be negative.
The universe set can be set via the parameter anno_universe
. By default it is the total items annotated to the whole DAG.
term_sim(dag, terms, method = "Sim_kappa", control = list(anno_universe = ...))
Definitions of the Jaccard, Dice and overlap coeffcients are similar. The Jaccard coeffcient is:
Jaccard(a,b)=|A∩B||A∪B|
The Dice coeffcient is:
Dice(a,b)=2∗|A∩B||A|+|B|
The overlap coeffcient is:
Overlap(a,b)=|A∩B|min{|A|,|B|}
Dice and Jaccard coeffcients have a relation of:
Jaccard=Dice2−Dice
The universe set can be set via the parameter anno_universe
. By default it is the total items annotated to the whole DAG.
term_sim(dag, terms, method = "Sim_Jaccard", control = list(anno_universe = ...)) term_sim(dag, terms, method = "Sim_Dice", control = list(anno_universe = ...)) term_sim(dag, terms, method = "Sim_Overlap", control = list(anno_universe = ...))
## R version 4.3.1 (2023-06-16)
## Platform: x86_64-apple-darwin20 (64-bit)
## Running under: macOS Ventura 13.2.1
##
## Matrix products: default
## BLAS: /Library/Frameworks/R.framework/Versions/4.3-x86_64/Resources/lib/libRblas.0.dylib
## LAPACK: /Library/Frameworks/R.framework/Versions/4.3-x86_64/Resources/lib/libRlapack.dylib; LAPACK version 3.11.0
##
## locale:
## [1] C/UTF-8/C/C/C/C
##
## time zone: Europe/Berlin
## tzcode source: internal
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] simona_1.1.3 knitr_1.44
##
## loaded via a namespace (and not attached):
## [1] sass_0.4.7 xml2_1.3.5 shape_1.4.6
## [4] stringi_1.7.12 digest_0.6.33 magrittr_2.0.3
## [7] evaluate_0.22 grid_4.3.1 RColorBrewer_1.1-3
## [10] iterators_1.0.14 circlize_0.4.15 fastmap_1.1.1
## [13] foreach_1.5.2 doParallel_1.0.17 rprojroot_2.0.3
## [16] jsonlite_1.8.7 GlobalOptions_0.1.2 promises_1.2.1
## [19] ComplexHeatmap_2.16.0 purrr_1.0.2 codetools_0.2-19
## [22] textshaping_0.3.7 jquerylib_0.1.4 shiny_1.6.0
## [25] cli_3.6.1 rlang_1.1.1 crayon_1.5.2
## [28] scatterplot3d_0.3-44 ellipsis_0.3.2 cachem_1.0.8
## [31] yaml_2.3.7 tools_4.3.1 parallel_4.3.1
## [34] memoise_2.0.1 colorspace_2.1-0 httpuv_1.6.11
## [37] GetoptLong_1.0.5 BiocGenerics_0.46.0 mime_0.12
## [40] vctrs_0.6.4 R6_2.5.1 png_0.1-8
## [43] matrixStats_1.0.0 stats4_4.3.1 lifecycle_1.0.3
## [46] stringr_1.5.0 S4Vectors_0.38.2 fs_1.6.3
## [49] IRanges_2.34.1 clue_0.3-65 cluster_2.1.4
## [52] ragg_1.2.6 pkgconfig_2.0.3 desc_1.4.2
## [55] later_1.3.1 pkgdown_2.0.7 bslib_0.5.1
## [58] Rcpp_1.0.11 glue_1.6.2 systemfonts_1.0.5
## [61] xfun_0.40 xtable_1.8-4 rjson_0.2.21
## [64] igraph_1.5.1 htmltools_0.5.6.1 rmarkdown_2.25
## [67] Polychrome_1.5.1 compiler_4.3.1