% \VignetteIndexEntry{An introduction to GOSemSim} % \VignetteDepends{org.Hs.eg.db,GO.db} % \VignetteSuggests{cluster} % \VignetteKeywords{GO Semantic Similarity Measurement} % \VignettePackage{GOSemSim} \documentclass[]{article} \usepackage{times} \usepackage{natbib} \usepackage{hyperref} \newcommand{\Rfunction}[1]{{\texttt{#1}}} \newcommand{\Robject}[1]{{\texttt{#1}}} \newcommand{\Rpackage}[1]{{\textit{#1}}} \newcommand{\Rmethod}[1]{{\texttt{#1}}} \newcommand{\Rfunarg}[1]{{\texttt{#1}}} \newcommand{\Rclass}[1]{{\textit{#1}}} \newcommand{\Rcode}[1]{{\texttt{#1}}} \newcommand{\R}{\textsf{R}} \newcommand{\GOSemSim}{\Rpackage{GOSemSim}} \newcommand{\Params}{\Rclass{Params}} \newcommand{\GOSet}{\Rclass{GOSet}} \newcommand{\GeneSet}{\Rclass{GeneSet}} \newcommand{\GeneClusterSet}{\Rclass{GeneClusterSet}} \newcommand{\term}[1]{\emph{#1}} \newcommand{\mref}[2]{\htmladdnormallinkfoot{#2}{#1}} \title{GO-terms Semantic Similarity Measures} \author{Guangchuang Yu \\ \\ Jinan University, Guangzhou, China} \begin{document} \bibliographystyle{plainnat} \maketitle <>= options(width=60) @ <>= ## w. verbose=TRUE library(GOSemSim) library(org.Hs.eg.db) library(GO.db) @ \section{Introduction} \indent Functional similarity of gene products can be estimated by controlled biological vocabularies, such as Gene Ontology (GO). GO comprises of three orthogonal ontologies, i.e. molecular function (MF), biological process (BP), and cellular component (CC). \\ \indent Four methods have been presented to determine the semantic similarity of two GO terms based on the annotation statistics of their common ancestor terms (Resnik\citep{philip_semantic_1999}, Jiang\citep{jiang_semantic_1997}, Lin\citep{lin_information-theoretic_1998} and Schlicker\citep{schlicker_new_2006}). Wang \citep{wang_new_2007} proposed a new method to measure the similarity based on the graph structure of GO. Each of these methods has its own advantages and weaknesses. The \textbf{GOSemSim} package \citep{g.yu_2010} is developed to compute semantic similarity among GO terms, sets of GO terms, gene products, and gene clusters, providing both five methods mentioned above. \\ \par \section{Semantic Similarity Measurement Based on GO} \indent The \textbf{GOSemSim} package contains functions to estimate semantic similarity of GO terms based on Resnik's, Lin's, Jiang and Conrath's, Rel's and Wang's method. Details about Resnik's, Lin's, and Jiang and Conrath's methods can be seen in \citep{lord_semantic_2003}, details about Rel's method can be seen in \citep{schlicker_new_2006}, and details about Wang's method can be seen in \citep{wang_new_2007}. \\ \indent Formally, a GO term A can be represented as $DAG_{A}=(A,T_{A},E_{A})$ where $T_{A}$ is the set of GO terms in $DAG_{A}$, including term A and all of its ancestor terms in the GO graph, and $E_{A}$ is the set of edges connecting the GO terms in $DAG_{A}$. \\ \indent To encode the semantics of a GO term in a measurable format to enable a quantitative comparison between two term's semantics, we firstly define the semantic value of term A as the aggregate contribution of all terms in $DAG_{A}$ to the semantics of term A, terms closer to term A in $DAG_{A}$ contribute more to its semantics. Thus, define the contribution of a GO term \textit{t} to the semantics of GO term A as the S-value of GO term \textit{t} related to term A. For any of term \textit{t} in $DAG_{A}=(A,T_{A},E_{A})$, its S-value related to term A. $S_{A}(\textit{t})$ is defined as: \begin{center} \[\left\{ \begin{array}{l} S_{A}(A)=1 \\ S_{A}(\textit{t})=\max\{w_{e} \times S_{A}(\textit{t}') | \textit{t}' \in childrenof(\textit{t}) \}$ if $\textit{t} \ne A \end{array} \right.\] \end{center} where $w_{e}$ is the semantic contribution factor for edge $e \in E_{A}$ linking term \textit{t} with its child term \textit{t}'. We defined term A contributes to its own as one. After obtaining the S-values for all terms in $DAG_{A}$, the semantic value of GO term A, SV(A), is calculated as: \begin{center} $SV(A)=\displaystyle\sum_{t \in T_{A}} S_{A}(t)$ \end{center} Given two GO terms A and B, the semantic similarity between these two terms, $GO_{A,B}$, is defined as: \begin{center} $S_{GO}(A,B)=\displaystyle\sum_{t \in T_{A} \cap T_{B}} \frac{S_{A}(t) + S_{B}(t)}{SV(A) + SV(B)}$ \end{center} where $S_{A}(\textit{t})$ is the S-value of GO term \textit{t} related to term A and $S_{B}(\textit{t})$ is the S-value of GO term \textit{t} related to term B. \\ \indent The method descriped above is proposed by Wang\citep{wang_new_2007}. The Wang's method determines the semantic similarity of two GO terms based on both the locations of these terms in the GO graph and their relations with their ancestor terms. \\ \indent The other four methods proposed by Resnik\citep{philip_semantic_1999}, Jiang\citep{jiang_semantic_1997}, Lin\citep{lin_information-theoretic_1998} and Schlicker\citep{schlicker_new_2006} are information content (IC) based, which depend on the frequencies of two GO terms involved and that of their closest common ancestor term in a specific corpus of GO annotations. Information content is defined as frequency of each term occurs in the corpus. At present, \textbf{GOSemSim} supports analysis on many species. We used Bioconductor package org.At.tair.db, org.Ag.eg.db, org.Bt.eg.db, org.Cf.eg.db, org.Gg.eg.db, org.Pt.eg.db, org.Sco.eg.db, org.EcK12.eg.db, org.EcSakai.eg.db, org.Dm.eg.db, org.Hs.eg.db, org.Pf.plasmo.db, org.Mm.eg.db, org.Ss.eg.db, org.Rn.eg.db, org.Mmu.eg.db, org.Ce.eg.db, org.Xl.eg.db, org.Sc.sgd.db and org.Dr.eg.db to calculate the information content of Arabidopsis, Anopheles, Bovine, Canine, Chicken, Chimp, Coelicolor, E coli strain K12 and strain Sakai, Fly, Human, Malaria, Mouse, Pig, Rat, Rhesus, Worm, Xenopus, Yeast and Zebrafish respectively. The information content will update regularly. \\ \indent As GO allow multiple parents for each concept, two terms can share parents by multiple paths. We take the mininmum p(t), where there is more than on shared parents. The \textit{$p_{ms}$} is defined as: \indent \begin{center} \textit{$p_{ms}(t1,t2)$} $= \displaystyle\min_{t \in S(t1,t2)} \{p(t)\})$ \end{center} Where S(t1,t2) is the set of parent terms shared by t1 and t2. \indent The similarity of Resnik's method is defined as: \begin{center} $sim(t1,t2) = -\ln p_{ms}(t1,t2)$ \end{center} \indent The similarity of Lin's is defined as: \begin{center} $sim(t1,t2)=\displaystyle\frac{2 \times \ln (p_{ms}(t1,t2))}{\ln p(t1) + \ln p(c2)}$ \end{center} \indent The similarity of Schlicker's method combine Resnik's and Lin's method is defined as: \begin{center} $sim(t1,t2)=\displaystyle\frac{2 \times \ln p_{ms}(t1,t2)}{\ln p(t1) + \ln p(p2)} \times (1-p_{ms}(t1,t2))$ \end{center} \indent The Jiang and Conrath's method define a semantic similarity as: \begin{center} $sim(t1,t2) = 1-\min(1, d(t1,t2))$ \end{center} where \begin{center} $d(t1,t2)= \ln p(t1) + \ln p(p2) - 2 \times \ln p_{ms}(t1,t2)$ \end{center} \indent In \textbf{GOSemSim}, on the basis of semantic similarity between GO terms, we can also compute semantic similarity among sets of GO terms, gene products, and gene clusters. \indent We implemented four methods which called \textit{max}, \textit{average}, \textit{rcmax}, and \textit{rcmax.avg} to combine semantic similarity scores of multiple GO terms. The similarities among gene products and gene clusters which annotated by multiple GO terms were also calculated by the same combine methods mentioned above. Given two GO terms sets $GO_{1}=\{go_{11},go_{12} \cdots go_{1m}\}$ and $GO_{2}=\{go_{21},go_{22} \cdots go_{2n}\}$, method \textit{max} calculate the maximum semantic similarity score over all pairs of GO terms between these two sets, method \textit{average} calcuate the average semantic similarity score over all pairs of GO terms. Similarities between GO terms form a matrix, and method \textit{rcmax} use the maximum of RowScore and ColumnScore as the similarity, where RowScore (or ColumnScore) is the average of maximum similarities on each row (or column). And method \textit{rcmax.avg} calculate the average of all maximum similarities on each row and column, and defined as: \begin{center} $Sim(GO1, GO2) = \frac{\displaystyle\sum_{1 \le i \le m} \max(Sim(\textit(go_{1i}), \textit(GO_{2}))) + \displaystyle\sum_{1 \le j \le n} \max(Sim(\textit(go_{2j}), \textit(GO_{1})))} {m+n}$ \end{center} \par \section{Functions} \indent A \Params{} stores a set of parameters for measuring semantic similarity. \Params{} containing parameters are ontology, organism, method, combine, and dropCodes. Parameter ontology specify which ontology were used in measurement, organism specifiy which GO Map were loaded for mapping Gene IDs to GO terms, dropCodes restrict evident codes when mapping Gene IDs to GO Terms, method specify which method to be used to measure the similarity and combine sepcify which combine method was used to combining semantic similarity scores. <>= params <- new("Params", ontology="MF", organism="human", method="Wang") @ A \GOSet{} stores two set of GO IDs. <>= go1 <- c("GO:0004022", "GO:0004024", "GO:0004023") go2 <- c("GO:0009055", "GO:0020037") gos <- new("GOSet", GOSet1=go1, GOSet2=go2) @ A \GeneSet{} containing two set of Gene IDs. <>= gs1 <- c("835", "5261","241", "994", "514", "533") gs2 <- c("578","582", "400", "409", "411") gs <- new("GeneSet", GeneSet1=gs1, GeneSet2=gs2) @ A \GeneClusterSet{} containing a list of gene clusters. <>= x <- org.Hs.egGO hsEG <- mappedkeys(x) clusters <- list(a=sample(hsEG, 20), b=sample(hsEG, 20), c=sample(hsEG, 20)) geneClusters <- new("GeneClusterSet", GeneClusters=clusters) @ Function \textit{sim} was designed to measuring semantic similarity among \GOSet{}, \GeneSet{} and \GeneClusterSet{}. <>= sim(gos,params) setCombineMethod(params)<-"rcmax.avg" sim(gos,params) sim(gs, params) sim(geneClusters, params) @ The old function calls which maybe more user friendly were also supported. They are wrapper functions of function \textit{sim} mentioning above. <>= goSim("GO:0004022", "GO:0005515", ont="MF", measure="Wang") go1 = c("GO:0004022","GO:0004024","GO:0004174") go2 = c("GO:0009055","GO:0005515") mgoSim(go1, go2, ont="MF", measure="Wang", combine="rcmax.avg") geneSim("241", "251", ont="MF", organism="human", measure="Wang", combine="rcmax.avg") mgeneSim(genes=c("835", "5261","241", "994"), ont="MF", organism="human", measure="Wang") clusterSim(gs1, gs2, ont="MF", organism="human", measure="Wang", combine="rcmax.avg") mclusterSim(clusters, ont="MF", organism="human", measure="Wang", combine="rcmax.avg") @ \section*{Session Information} The version number of R and packages loaded for generating the vignette were: \begin{verbatim} <>= sessionInfo() @ \end{verbatim} \bibliography{GOSemSim} \end{document}