% -*- mode: noweb; noweb-default-code-mode: R-mode; -*- %\VignetteIndexEntry{How to use MiRaGE Package} %\VignetteKeywords{miRNA, gene expression} %\VignetteDepends{MiRaGE} %\VignettePackage{MiRaGE} %documentclass[12pt, a4paper]{article} \documentclass[12pt]{article} \usepackage{amsmath} \usepackage{hyperref} %\usepackage[authoryear,round]{natbib} \textwidth=6.2in \textheight=8.5in %\parskip=.3cm \oddsidemargin=.1in \evensidemargin=.1in \headheight=-.3in \newcommand{\scscst}{\scriptscriptstyle} \newcommand{\scst}{\scriptstyle} \newcommand{\Rfunction}[1]{{\texttt{#1}}} \newcommand{\Robject}[1]{{\texttt{#1}}} \newcommand{\Rpackage}[1]{{\textit{#1}}} \author{Y-h. Taguchi} \begin{document} \title{How to use MiRaGE Package} \maketitle \tableofcontents \section{Introduction} This document describes briefly how to use \verb+MiRaGE+ package in order to infer the target gene regulation by miRNA, based upon target gene expression. MiRaGE is based upon the algorithm proposed in \cite{LNCS,LNBI,IJMS}. Basically, its function is the same as the MiRaGE Server. In order to infer the target gene regulation by miRNAs, we made use of target gene (mRNA) expression. Suppose $x_{gs}$ is the expression of the $g$th gene in the $s$th sample. Then $P$-value to measure the amount of target gene regulation by the $m$th miRNAs is computed by several statistical test. More detailed and comprehensive explanations can be found in \cite{iConceptpress}. \section{Background} miRNA is short non-coding RNA (ncRNA) which is believed to degrade target genes. Target genes are believed to be decided by seed match between 7-mer at 5' untranslated region (UTR) of miRNA and 3' UTR of target mRNAs. However because of huge number of miRNAs (c.a. 1000) and the huge number of target genes (c.a. hundreds) of each miRNA, it is not easy to experimentally decide which miRNA regulates target genes. MiRaGE infers target gene regulation from target gene expression and computationally predicted target gene table. It gives the rejection probability to reject null hypothesis that target genes of a specific miRNA are equally regulated as other genes. When $t$-test is employed, $P$-value is $$ P(S^{ss'}_m > S'^{ss'}_m) $$ or $$ P(S^{ss'}_m < S'^{ss'}_m) $$ where $P$ is the rejection probability of null hypothesis $S_m= S'_m$ when the alternative hypothesis is either $S_m > S'_m$ or $S_m < S'm$. $S_m$ and $S'_m$ are the test variable to measure the target gene regulation by the $m$th miRNA. When $P$-values are computed via $t$-test, $S_m$ is the mean gene expression logarithmic ratio of the $m$th miRNA's target genes, i.e., $$ S^{ss'}_m = \frac{1}{N(G_m)}\sum_{g \in G_m} \log \frac{x_{gs}}{x_{gs'}} $$ where $G_m$ is the set of the $m$th miRNA's target genes and $N(G_m)$ is the total number of genes in $G_m$. $S'_m$ is the mean expression logarithmic ratio of genes not targeted by the $m$th miRNA but any other miRNAs and is defined as $$ S'^{ss'}_m = \frac{1}{N(G'_m)}\sum_{g \in G'_m} \log \frac{x_{gs}}{x_{gs'}} $$ where $G'_m$ is the set of the $m$th miRNA's target genes and $N(G'_m)$ is the total number of genes in $G'_m$. On the other hands, if $P$-values are computed by Wilcoxon rank sum test, they are $$ P\left(U_m^{ss'}>\frac{N(G_m)N(G'_m)}{2}\right) $$ or $$ P\left(U_m^{ss'}<\frac{N(G_m)N(G'_m)}{2}\right) $$ which are the rejection probabilities of null hypothesis $U= \frac{N(G_m)N(G'_m)}{2}$ when the alternative hypothesis is either $U> \frac{N(G_m)N(G'_m)}{2}$ or $U < \frac{N(G_m)N(G'_m)}{2}$. Here $U^{ss'}$ is the test variable, $$ U_m^{ss'} = R^{ss'}_m - \frac{N(G_m) (N(G_m)+1)}{2} $$ Here $$ R^{ss'}_m = \sum_{g \in G_m} R^{ss'}\left (\log \frac{x_{gs}}{x_{gs'}} \right) $$ where $R^{ss'}( \ldots)$ is the rank order of the $m$th genes logarithmic ratio among all of considered genes. Alternatively, Kolmogorov-Smirnov test can be employed. In this case, test variable is $$ D_m^{ss'} = \sup_g \left(F_m(x_{gs}) - F'_m(x_{gs})\right) $$ or $$ {D'}_m^{ss'} = \sup_g \left(F'_m(x_{gs}) - F_m(x_{gs})\right) $$ where $$ F^{ss'}_m(x_{gs}) = \frac{1}{N(G_m)} \sum_{g' \in G_m} \Theta \left(\log \frac{x_{gs}}{x_{gs'}} - \log \frac{x_{g's}}{x_{g's'}}\right) $$ and $$ {F'}^{ss'}_m(x_{gs}) = \frac{1}{N(G'_m)} \sum_{g' \in G'_m} \Theta \left(\log \frac{x_{gs}}{x_{gs'}} - \log \frac{x_{g's}}{x_{g's'}}\right) $$ where $\Theta(x)$ is the step function, $$ \Theta(x) = \left \{\begin{array}{ll} 1 & x>0 \\ 0 & x<0 \end{array} \right. $$ Then $P$-value is computed via $P(D_m^{ss'}>0)$ or $P({D'}_m^{ss'}>0)$ under the null hypothesis $D_m^{ss'}=0$ or ${D'}_m^{ss'}=0$. Since the target genes table is generated by the simple seed match, MiRaGE does not need any other external programs to obtain target gene table. Another advantage is the exclusion of mRNA targeted by no miRNAs. This enables us more accurate prediction of target gene regulation by miRNAs. \section{Quick start} %\begin{Sinput} %R> library(MiRaGE) ##load the MiRaGE package %R> data(x_gena) ##load sample data %R> result <- MiRaGE(x_gene) ## excute MiRaGE function %\end{Sinput} %<>= <<>>= library(MiRaGE) data(gene_exp) library(Biobase) result <- MiRaGE(gene_exp,species="HS") @ Then \verb+result$P0+ and \verb+result$P1+ include $P$-values for upregulation and downregulation of target genes by miRNAs, respectively. The definition ``up" or ``down" depends upon the order of columns of expression data in \verb+gene_exp+ (see below). \underline{\bf Caution} I strongly recommend user to use \verb+location+="web" option, since it will be most frequently updated. Default setting requires experimental package \Rpackage{miRNATarget} (see Sec. \ref{miRNATarget}). Data set on the web can be stored for the later usage, too (see below). <>= result <- MiRaGE(gene_exp,location="web",species="HS") @ \section{Data Structure} \subsection{Input: target gene expression} In order to execute analysis, you need ExpressionSet objects which stores target gene expression in it. In order to see this, it is easier to see sample data \verb+gene_exp+ as follows. <<>>= data(gene_exp) gene_exp @ The above example displays the ExpressionSet object, \verb+gene_exp+. As you can see, each row corresponds to each gene and each column corresponds to each sample (experiment). \verb+gene_id+ in featureData must includes gene id. \verb+sample_name+ in phenoData must include sample names. Since \verb+MiRaGE+ package tries to compare two distinct states, you need at least a set of gene expression corresponding to each of them. In \verb+gene_exp+ data, we have two biological replicates of negative control and the results one day after treatment. Thus, the 1st and 2nd columns of expression data are named as \verb+neg.1+ and \verb+neg.2+, respectively (This means ``negative control 1" and ``negative control 2", respectively). The 3rd and 4th columns of expression data corresponds to the two biological replicates one day after the treatment. Thus, they are names as \verb+day1.1+ and \verb+day1.2+, respectively. These column names which express distinct samples keep some flexibilities but must have the form $ group.n$, where $group$ corresponds to either of sample groups and $n$ must be integer starting from 1. This means, it you have $N$ biological replicates for the first group (typically it includes un-treated or negative control samples) names as $groupA$ and $M$ for the second group (typically it includes treated samples) names as $groupB$, data structure of ExpressioSet which stores target gene expression is, \begin{verbatim} sample_name: groupA.1 groupA.1 ... groupA.N groupB.1 groupB.2 ... groupB.M gene_id : gene1,gene2,gene3.... \end{verbatim} \verb+gene_id+ is much easier. It can includes any of gene id which can be treated by \Rfunction{MiRaGE}. They can be a mixture of the different types of gene ids. In this case, only gene expression having gene id specified when \Rfunction{MiRaGE} is called are treated as target genes. The easiest way to generate ExpressionSet which include target gene expression may be importing files including gene expression using standard \verb+R+ function, <<>>= x_gene <- read.csv(system.file("extdata/x_all_7a.csv",package="MiRaGE"),sep="\t") x_gene[101:103,] @ As can be seen, the first column includes gene id, which is "refseq" here, and the second to the fifth columns include gene expression. Data frame \verb+x_gene+ can be transformed to ExpressionSet objects \verb+gene_exp+ as <<>>= gene_exp <- new("ExpressionSet",expr=data.matrix(x_gene[,-1])) fData(gene_exp)[["gene_id"]] <- x_gene[,1] pData(gene_exp)[["sample_name"]] <- colnames(x_gene)[-1] @ For users' convenience, we have places a file \verb+x_all_7a.csv+ under \verb+csv+ directory. Please refer to this file for the preparation of files including target gene expression. \subsection{Output: $P$-values} As mentioned in the above, output of \Rfunction{MiRaGE} is a list which includes two dataframes named as \verb+P0+ and \verb+P1+, respectively. \verb+P0+ includes the rejection probabilities that the target gene expression in the first sample group is less than that in the second group. This means, smaller $P$-values indicate the target gene expression in the first sample group is more likely less than the second sample groups. Inversely, \verb+P1+ includes the rejection probabilities that the target gene expression in the second sample is less than that in the first group. Thus, smaller $P$-values indicate target gene expression in the second group is more likely less than the first groups. <<>>= result$P1[1:3,] @ In the above, we have shown the first three lines in the dataframe \verb+result$P1+. Since these are small, target genes of these three miRNAs is possibly expressive in the second group. In the first column of \verb+result$P0+ and \verb+result$P1+, names of considered miRNAs are listed. The number of miRNAs considered varys dependent upon the argument \verb+conv+ of \Rfunction{MiRaGE}. The second column includes $P$-values attributed to each miRNA. Dependent upon argument \verb+method+, the number of columns which store $P$-values may change (see below). \section{Example} \subsection{Example1: non-differentiated vs differentiated ES cell} In this section, we demonstrate how to infer target gene regulation via \Rfunction{MiRaGE}. First we import data set from experimental package \Rpackage{humanStemCell} <<>>= require(humanStemCell) data(fhesc) @ In this data set, human stem cells were assayed using Affymetrix 133plus 2 arrays. There were six arrays, three were biological replicates for undifferentiated cells, the other three were biological replicates for differentiated cells. In order to analyze this set, we modify ExpressionSet \verb+fhesc+ as <<>>= pData(fhesc)[["sample_name"]] <- c("neg.1","neg.2","neg.3", "pos.1","pos.2","pos.3") fData(fhesc)[["gene_id"]] <-rownames(exprs(fhesc)) @ Then, first three are designated as non-differentiated ES cell and the later three are differentiated ES cell. Obtaining P-values is easy, <<>>= require(MiRaGE) result <- MiRaGE(fhesc,species="HS",ID="affy_hg_u133a_2") @ Using the results, we can list miRNAs whose target genes are upregulated in the later (i.e., differentiated ES cell) group with P-values. <<>>= result$P0[order(result$P0[,2])[1:5],] @ Since miRNAs are believed to suppress target genes, these miRNAs are supposed to be upregulated in the former (i.e., non-differentiated ES cell) group. \subsection{Example 2: Universal Human Reference RNA vs brain} In this section, we demonstrate how to infer target gene regulation via \Rfunction{MiRaGE} using another example. First we import data set from experimental package \Rpackage{beadarrayExampleData} <<>>= require(beadarrayExampleData) data(exampleBLData) data(exampleSummaryData) @ The data in this package are a subset of the MAQC bead-level data available in the beadarrayUseCases package. Bead-level refers to the availability of intensity and location information for each bead on each BeadArray in an experiment. In this dataset, BeadArrays were hybridized with either Universal Human Reference RNA (UHRR, Stratagene) or Brain Reference RNA (Ambion) as used inthe MAQC project. This object is a representation of the bead-level data for 2 arrays and was created by the beadarray package. Since this is two color array, and the number of columns of expression must be the number of columns of expression data MUST be the length of \verb+sample_name+, we omit later half of samples and employ only the first twelve samples, for simplicity. <<>>= vv <- exampleSummaryData[,1:12] fData(vv)[["gene_id"]] <- fData(exampleSummaryData)[["IlluminaID"]] pData(vv)[["sample_name"]] <- c("neg.1","neg.2","neg.3","neg.4", "neg.5","neg.6","brain.1","brain.2","brain.3","brain.4","brain.5","brain.6") result <- MiRaGE(vv,species="HS",ID="illumina_humanwg_6_v3") @ Then we can list miRNAs whose target genes are upregulated in negative control, i.e., miRNAs which are expected to be upregulated in brain as follows. <<>>= result$P1[order(result$P1[,2])[1:5],] @ \section{Rapid use \& Off line use} Although the default value of \verb+location+ is "local", when \verb+location+="web", \Rfunction{MiRaGE} every time tries to access \verb+MiRaGE Server+\footnote{http://www.granular.com/DATA2/} to download target gene tables, gene id conversion table, and miRNA conservation table. It is a time consuming process. Especially, since the target gene table is huge, it may take a few minutes. It may not be often to use \Rfunction{MiRaGE} iteratively many times, we offer the method to avoid ``every time download". \subsection{Suppressing downloading} In \Rfunction{MiRaGE}, we offer the option to suppress downloading. If you repeatedly use \Rfunction{MiRaGE} with keeping either \verb+species+, \verb+ID+, or \verb+conv+ unchanged, you can suppress time consuming download process by specifying either \verb+species_force+, \verb+ID_force+, or \verb+conv_force+ as \verb+FALSE+ (Defaults for these are \verb+TRUE+). \underline{\bf Caution} Do not omit the arguments either \verb+species+, \verb+ID+, or \verb+conv+ if they differ from defaults, even if they are not modified during iterative usage and either \verb+species_force+, \verb+ID_force+, or \verb+conv_force+ is \verb+FALSE+. They are used for other purposes than specifying what should be downloaded. \subsection{Save \& load tables} More advanced and convenient way is to save the objects storing target gene tables, gene id conversion table, and miRNA conservation table. The names of objects are, \begin{itemize} \item{\verb+TBL2+} : Target gene tables \item{\verb+id_conv+} : Gene id conversion table \item{\verb+conv_id+} : MiRNA conservation table \end{itemize} Thus, for example \verb+TBL2+ is saved as <<>>= save(file="TBL2",TBL2) @ you can use it later by loading as <<>>= load("TBL2") @ Then you can execute \Rfunction{MiRaGE} with specifying \verb+species_force=F+ as <>= result <- MiRaGE(...,species_force=F) @ Now, you can skip time consuming download processes for the target gene table. Similar procedures are possible for \verb+id_conv+ and \verb+conv_id+, too. Execute \Rfunction{MiRaGE}, save downloaded tables, and use the tables later by loading them when these arguments take same values. \subsection{miRNATarget package} \label{miRNATarget} One can also install experimental package \Rpackage{miRNATarget} instead of the usage of web. Once you install experimental package \Rpackage{miRNATarget}, you will never be required to access to internet. <<>>= library(MiRaGE) data(gene_exp) library(Biobase) result <- MiRaGE(gene_exp,species="HS") @ \subsection{Generation of tables from scratch} I have also prepared functions which generate \verb+TBL2+, \verb+id_conv+ and \verb+conv_id+ from scratch. Usually, user do not need them since prepared tables can be obtained from the web or as experimental package as mentioned above. \verb+TBL2_HS+ can be saved in the current directly by executing <>= TBL2_HS_gen() @ and \verb+TBL2_MM+ can saved in the current directly by executing <>= TBL2_MM_gen() @ \verb+id_conv+ for mouse can saved in the current directly by executing <>= id_conv_gen(SP="MM") @ and \verb+id_conv+ for human can saved in the current directly by executing <>= id_conv_gen(SP="HS") @ \verb+HS_conv_id+ can saved in the current directly by executing <>= HS_conv_id() @ and \verb+MM_conv_id+ can saved in the current directly by executing <>= MM_conv_id() @ However, basically, execution of some of them are very time consuming. It is highly discouraged to build tables from scratch. It is much better to use prepared tables. \section{Multiple comparison correction} Obtained $P$-values are definitely underestimated, i.e., even if $P<0.05$, this does not mean the rejection probability is less than 0.05. If one prefers to use adjusted $P$-values, we recommend to use \Rfunction{p.adjust} with parameter of \verb+BH+, as <<>>= p.adjust(result$P1[,2],method="BH") @ Then we can see which $P$-values are really significant, e.g., less than $0.05$. In addition to this, it will allow us to evaluate which miRNAs really regulate target genes, e.g., <<>>= result$P1[,1][p.adjust(result$P1[,2],method="BH")<0.05] @ \Robject{x\_gene} is the transfection expreiments of let-7a, it is reasnable that only a few miRNAs including let-7a have significant $P$-values. \section*{Appendix:Arguments} Functionality of \Rfunction{MiRaGE} changes dependent upon the values of arguments. In this section, we will try to explain how the functionality of \Rfunction{MiRaGE} changes. \subsection*{species} This specifies target species. Considered miRNAs are based upon miRBase\footnote{http://www.mirbase.org}. Rel. 20. At the moment, supported species are human (``HS") and mouse (``MM"). \Rfunction{MiRaGE} downloads corresponding target gene table (named \verb+TBL2+) from MiRaGE Server. Default is ``MM". \subsection*{ID} This specifies gene ID. Default is ``refseq". If ID is not ``refseq", \Rfunction{MiRaGE} downloads corresponding gene id conversion table (called \verb+ID+) from RefSeq to specified gene ID from MiRaGE Server. Supported gene IDs are, \bigskip \begin{tabular}{rll} \hline \multicolumn{3}{c}{common} \\ & \verb+ID+ & description \\\hline 1 & ensembl\_gene\_id & Ensembl Gene ID \\ 2 & ensembl\_transcript\_id & Ensembl Transcript ID \\ 3 & ensembl\_peptide\_id & Ensembl Protein ID \\ 4 & ensembl\_exon\_id & Ensembl Exon ID \\ 5 & ccds & CCDS ID \\ 6 & embl & EMBL (Genbank) ID \\ 7 & entrezgene & EntrezGene ID \\ 8 & merops & MEROPS ID \\ 9 & pdb & PDB ID \\ 10 & protein\_id & Protein (Genbank) ID \\ 11 & refseq\_peptide & RefSeq Protein ID [e.g. NP\_001005353] \\ 12 & rfam & Rfam ID \\ 13 & rfam\_transcript\_name & Rfam transcript name \\ 14 & ucsc & UCSC ID \\ 15 & unigene & Unigene ID \\ 16 & uniprot\_sptrembl & UniProt/TrEMBL Accession \\ 17 & uniprot\_swissprot & UniProt/SwissProt ID \\ 18 & uniprot\_swissprot\_accession & UniProt/SwissProt Accession \\ 19 & uniprot\_genename & UniProt Gene Name \\ 20 & uniprot\_genename\_transcript\_name & Uniprot Genename Transcript Name \\ 21 & wikigene\_name & WikiGene Name \\ 22 & wikigene\_id & WikiGene ID \\ 23 & efg\_agilent\_sureprint\_g3\_ge\_8x60k & Agilent SurePrint G3 GE 8x60k probe \\ 24 & efg\_agilent\_wholegenome\_4x44k\_v1 & Agilent WholeGenome 4x44k v1 probe \\ 25 & efg\_agilent\_wholegenome\_4x44k\_v2 & Agilent WholeGenome 4x44k v2 probe \\ 26 & codelink & Codelink probe \\ 27 & phalanx\_onearray & Phalanx OneArray probe \\ 28 & smart & SMART ID \\ 29 & pfam & PFAM ID \\ 30 & tigrfam & TIGRFam ID \\ 31 & interpro & Interpro ID \\ \hline \end{tabular} \bigskip \begin{tabular}{rll} \hline \multicolumn{3}{c}{human}\\ & \verb+ID+ & description \\\hline1 & hgnc\_id & HGNC ID(s) \\ 2 & hgnc\_symbol & HGNC symbol \\ 3 & hgnc\_transcript\_name & HGNC transcript name \\ 4 & affy\_hc\_g110 & Affy HC G110 probeset \\ 5 & affy\_hg\_focus & Affy HG FOCUS probeset \\ 6 & affy\_hg\_u133\_plus\_2 & Affy HG U133-PLUS-2 probeset \\ 7 & affy\_hg\_u133a\_2 & Affy HG U133A\_2 probeset \\ 8 & affy\_hg\_u133a & Affy HG U133A probeset \\ 9 & affy\_hg\_u133b & Affy HG U133B probeset \\ 10 & affy\_hg\_u95av2 & Affy HG U95AV2 probeset \\ 11 & affy\_hg\_u95b & Affy HG U95B probeset \\ 12 & affy\_hg\_u95c & Affy HG U95C probeset \\ 13 & affy\_hg\_u95d & Affy HG U95D probeset \\ 14 & affy\_hg\_u95e & Affy HG U95E probeset \\ 15 & affy\_hg\_u95a & Affy HG U95A probeset \\ 16 & affy\_hugenefl & Affy HuGene FL probeset \\ 17 & affy\_huex\_1\_0\_st\_v2 & Affy HuEx 1\_0 st v2 probeset \\ 18 & affy\_hugene\_1\_0\_st\_v1 & Affy HuGene 1\_0 st v1 probeset \\ 19 & affy\_u133\_x3p & Affy U133 X3P probeset \\ 20 & agilent\_cgh\_44b & Agilent CGH 44b probe \\ 21 & illumina\_humanwg\_6\_v1 & Illumina HumanWG 6 v1 probe \\ 22 & illumina\_humanwg\_6\_v2 & Illumina HumanWG 6 v2 probe \\ 23 & illumina\_humanwg\_6\_v3 & Illumina HumanWG 6 v3 probe \\ 24 & illumina\_humanht\_12 & Illumina Human HT 12 probe \\ \hline \end{tabular} \bigskip \begin{tabular}{rll} \hline \multicolumn{3}{c}{mouse}\\ \verb+ID+ & description \\\hline 1 & fantom & Fantom ID \\ 2 & ipi & IPI ID \\ 3 & mgi\_id & MGI ID \\ 4 & mgi\_symbol & MGI symbol \\ 5 & mgi\_transcript\_name & MGI transcript name \\ 6 & affy\_mg\_u74a & Affy mg u74a probeset \\ 7 & affy\_mg\_u74av2 & Affy mg u74av2 probeset \\ 8 & affy\_mg\_u74b & Affy mg u74b probeset \\ 9 & affy\_mg\_u74bv2 & Affy mg u74bv2 probeset \\ 10 & affy\_mg\_u74c & Affy mg u74c probeset \\ 11 & affy\_mg\_u74cv2 & Affy mg u74cv2 probeset \\ 12 & affy\_moe430a & Affy moe430a probeset \\ 13 & affy\_moe430b & Affy moe430b probeset \\ 14 & affy\_moex\_1\_0\_st\_v1 & Affy MoEx probeset \\ 15 & affy\_mogene\_1\_0\_st\_v1 & Affy MoGene probeset \\ 16 & affy\_mouse430\_2 & Affy mouse430 2 probeset \\ 17 & affy\_mouse430a\_2 & Affy mouse430a 2 probeset \\ 18 & affy\_mu11ksuba & Affy mu11ksuba probeset \\ 19 & affy\_mu11ksubb & Affy mu11ksubb probeset \\ 20 & illumina\_mousewg\_6\_v1 & Illumina MouseWG 6 v1 probe \\ 21 & illumina\_mousewg\_6\_v2 & Illumina MouseWG 6 v2 probe \\ \hline \end{tabular} \bigskip Requirements for supporting any other gene IDs are welcomed. \subsection*{method} This specifies how to treat replicates. if \verb+method+ is ``mean", then averaged gene expression is attributed to each gene. If it is ``mixed", they are used for statistical test as it is. This means, the number of target genes attributed to each miRNAs is as many as the number of replicates. If "one\_by\_one" is specified, all of combinations between the two groups, i.e., \begin{center} groupA.1 $\times$ groupB.1, groupA.1 $\times$ groupB.2, \ldots , groupA.N $\times$ groupB.M. \end{center} are condiered. Thus, in this case, both \verb+P0+ and \verb+P1+ have $1 + N \times M$ columns, the later $N \times M$ includes $P$-value for each of combinations. Default is ``mean". \subsection*{test} This specifies the statititical methods to evaluate siginificance of reglation of target genes. Supported are ``ks" (Kolmogorov-Smirnov test), ``t" (t-test), and "wilcox" (Wilcoxon test). These are performed by standard \verb+R+ functions, \Rfunction{ks.test}, \Rfunction{t.test}, and \Rfunction{wilcox.test}, respectively. Default is "ks". \subsection*{conv} This specifies how well considered miRNAs must be conserved. Supported are ``conserved", ``weak\_conserv" and ``all". Baed upon TargetScan 7.2 \footnote{http://www.targetscan.org}, they correspond to broadly conserved, conserved, and others. For more detail, plese colusult with TargetScan. Default is ``conserved". \subsection*{Force download or not} \verb+species_force+, \verb+ID_force+, and \verb+conv_force+ spefify if target gene table, gene id conversion table, and miRNA conservation table are forced to be donwloaded. Dafult is \verb+T+. If some of them have already been downloded and one would like to use it as it is, please specify they are \verb+F+. \begin{thebibliography}{99} \bibitem{LNCS} Y-h Taguchi, Jun Yasuda, 2010, Inference of Gene Expression Regulation via microRNA Transfection, ICIC2010, Proceedings, Springer, 6215, 672-679. \bibitem{LNBI} Y-h Taguchi, Jun Yasuda, 2012, MiRaGE: Inference of Gene Expression Regulation via MicroRNA Transfection II, ICIC2011, Proceedings, Springer, 6840,192-135 \bibitem{IJMS} M. Yoshizawa, Y-h. Taguchi, Jun Yasuda, 2011, Inference of Gene Regulation via miRNAs During ES Cell Differentiation Using MiRaGE Method, Int. J. Mol. Sci., 12[12]:9265-9276 \bibitem{iConceptpress} Taguchi, Y-h. (2013). Inference of Target Gene Regulation by miRNA via MiRaGE Server. Introduction to Genetics - DNA Methylation, Histone Modification and Gene Regulation. ISBN: 978-1477554-94-4. iConcept Press. Retrieved from http://www.iconceptpress.com/books/ IntroductionToGeneticsDNAMethylationHistoneModificationAndGeneRegulation/ \end{thebibliography} \end{document}