twoGOSim: Protein Similarity Calculation based on Gene Ontology (GO)...

View source: R/par-02-parGOSim.R

twoGOSimR Documentation

Protein Similarity Calculation based on Gene Ontology (GO) Similarity

Description

This function calculates the Gene Ontology (GO) similarity between two groups of GO terms or two Entrez gene IDs.

Usage

twoGOSim(
  id1,
  id2,
  type = c("go", "gene"),
  ont = c("MF", "BP", "CC"),
  organism = "human",
  measure = "Resnik",
  combine = "BMA"
)

Arguments

id1

Character vector. When length > 1: each element is a GO term; when length = 1: the Entrez Gene ID.

id2

Character vector. When length > 1: each element is a GO term; when length = 1: the Entrez Gene ID.

type

Input type of id1 and id2, 'go' for GO Terms, "gene" for gene ID.

ont

Default is "MF", can be one of "MF", "BP", or "CC" subontologies.

organism

Organism name. Default is "human", can be one of "anopheles", "arabidopsis", "bovine", "canine", "chicken", "chimp", "coelicolor", "ecolik12", "ecsakai", "fly", "human", "malaria", "mouse", "pig", "rat", "rhesus", "worm", "xenopus", "yeast" or "zebrafish". Before specifying the organism, please install the corresponding genome wide annotation data package for the selected organism.

measure

Default is "Resnik", can be one of "Resnik", "Lin", "Rel", "Jiang" or "Wang".

combine

Default is "BMA", can be one of "max", "average", "rcmax" or "BMA" for combining semantic similarity scores of multiple GO terms associated with proteins.

Value

Similarity value.

Author(s)

Nan Xiao <https://nanx.me>

See Also

See parGOSim for protein similarity calculation based on Gene Ontology (GO) semantic similarity. See parSeqSim for paralleled protein similarity calculation based on Smith-Waterman local alignment.

Examples

## Not run: 

# Be careful when testing this since it involves GO similarity computation
# and might produce unpredictable results in some environments

library("GOSemSim")
library("org.Hs.eg.db")

# By GO terms
go1 <- c("GO:0004022", "GO:0004024", "GO:0004023")
go2 <- c("GO:0009055", "GO:0020037")
twoGOSim(go1, go2, type = "go", ont = "MF", measure = "Wang")

# By Entrez gene id
gene1 <- "1956" # EGFR
gene2 <- "2261" # FGFR3
twoGOSim(gene1, gene2, type = "gene", ont = "BP", measure = "Lin")

## End(Not run)

protr documentation built on Nov. 2, 2023, 6:04 p.m.