Title: | Identification and Analysis of ceRNA Regulation |
---|---|
Description: | Provides several functions to identify and analyse miRNA sponge, including popular methods for identifying miRNA sponge interactions, two types of global ceRNA regulation prediction methods and four types of context-specific prediction methods( Li Y et al.(2017) <doi:10.1093/bib/bbx137>), which are based on miRNA-messenger RNA regulation alone, or by integrating heterogeneous data, respectively. In addition, For predictive ceRNA relationship pairs, this package provides several downstream analysis algorithms, including regulatory network analysis and functional annotation analysis, as well as survival prognosis analysis based on expression of ceRNA ternary pair. |
Authors: | Mengying Zhang,Yongsheng Li,Juan Xu*,Xia Li* |
Maintainer: | Mengying Zhang <[email protected]> |
License: | GPL-3 |
Version: | 2.1.3 |
Built: | 2024-11-12 06:30:42 UTC |
Source: | CRAN |
identifying miRNA sponge interactions .
CeRNASeek considers five method to identify ceRNA crosstalk, including two types of global ceRNA regulation prediction methods(ratio based, we termed ratio and hypergeometric test based, termed HyperT) and three types of context-specific prediction methods (hypergeometric test plus coexpression, termed HyperC, sensitivity correlation-based method (SC)and conditional mutual information (CMI)-based methods and cernia method). Ratio-based prediction: this method ranked the candidate genes by the proportion of common miRNAs,For instance, we need to identify the ceRNA partners for gene i from all candidate gene sets S, and the ratio is calculated as the intersection of miRNAi and miRNAj divided by miRNAj Where miRNAi is the miRNA set that regulated gene i and miRNAj is the miRNA set that regulated gene j.
Hypergeometric test-based prediction-HyperT: it is usually used the hypergeometric to evaluate whether two genes were coregulated by miRNAs,This statistic test computed the significance of common miRNAs for each RNA pairs.
Hypergeometric test combined with coexpression-based prediction-HyperC: To discover the active ceRNA-ceRNA regulatory pairs in a specific context, the commonly used method is to using the coexpression principle to filter the ceRNA-ceRNA regulation identified based on the above two global methods,this method integrated context-specific gene expression profile data sets. The Pearson correlation coefficient (R) of each candidate ceRNA regulatory pairs identified was calculated.
SC-based prediction: Another method introduced the expression profile data of shared miRNAs, and uses partial correlation to calculate ceRNA interaction pairs. Then, the Sensitivity correlation(SC) of miRNA-M, for the corresponding candidate ceRNA pair is calculated.
CMI-based methods: CMI is widely used to identify the RNA-RNA correlations, given the value of miRNAs. Hermes is a widely used method, which predicts ceRNA interactions from expression profiles of candidate RNAs and their common miRNA regulators using CMI, as in Sumazin et al.Firstly computed the significance of difference of I(miRNAi;T2|T1) and I(miRNAi;T2),where miRNAi represents the ith miRNA in the miRNA shared by the two target genes. Then Shuffled condition target's expression in 1000 times. The final significance based on Fisher’s method was calculated.For more information, please refer to the article by Sumazin P et al.
cernia method cernia method is implemented based on the following seven scores: 1. The fraction of common miRNAs; 2. The density of the MREs for all shared miRNAs; 3. The distribution of MREs of the putative ceRNAs; 4. The relation between the overall number of MREs for a putative ceRNA, compared with the number of miRNAs that yield these MREs; 5. The density of the hybridization energies related to MREs for all shared miRNAs; 6. The DT-Hybrid recommendation scores; and 7. The pairwise Pearson correlation between putative ceRNA expressions from selected tissue.
Xu J , Feng L , Han Z , et al. Extensive ceRNA–ceRNA interaction networks mediated by miRNAs regulate development in multiple rhesus tissues[J]. Nucleic Acids Research, 2016:gkw587.
Sumazin P , Yang X , Chiu H S , et al. An Extensive MicroRNA-Mediated Network of RNA-RNA Interactions Regulates Established Oncogenic Pathways in Glioblastoma[J]. Cell, 2011, 147(2):0-381.
Paci P , Colombo T , Farina L . Computational analysis identifies a sponge interaction network between long non-coding RNAs and messenger RNAs in human breast cancer[J]. BMC Systems Biology, 2014, 8(1):83.
Zhang Y , Xu Y , Feng L , et al. Comprehensive characterization of lncRNA-mRNA related ceRNA network across 12 major cancers[J]. Oncotarget, 2014, 7(39):64148-64167.
Sardina D S , Alaimo S , Ferro A , et al. A novel computational method for inferring competing endogenous interactions[J]. Briefings in Bioinformatics, 2016:bbw084.
The ceRNA function is used to identify the miRNAs that co-regulate the gene pairs of interest. The users can input the genes of interest. Otherwise, all the paired genes were calculated the number of co-regulating miRNAs.
ceRNA.app(miRtar, targetce = NULL)
ceRNA.app(miRtar, targetce = NULL)
miRtar |
A data frame representing the relationship between miRNA and target. The data frame contains the name of the miRNA and target regulatory relationship. |
targetce |
a character string (vector) specifying candidate target name to analyse (default (targetce = NULL)). |
Note:All the arguments without default value must be assigned.
A list of identified miRNAs that co-regulate the gene pairs containing following components:
ceRNA
identify the miRNAs that co-regulate the gene pairs of interest,List all possible
ceRNA interactions,a 4 columns dataframe as following:
targetce
represented target names,respectively.
anotherce
names of modulators that another possible target(modulators) constitutes a ceRNA interaction relation.
miRNAs
names of miRNA shared by two targets.
miRNAs_num
number of miRNAs shared by two targets.
miR_l
Number of miRNAs interacting with each target in the input file.
##Here we take the regulatory relationship between six genes and 71 miRNAs. ceRNA.app(dataset[["miRtar"]],targetce=NULL)
##Here we take the regulatory relationship between six genes and 71 miRNAs. ceRNA.app(dataset[["miRtar"]],targetce=NULL)
This function predicts ceRNA interactions by ratio or hypergeometric test.
ceRNA.basic(miRtar, targetce = NULL, method = "ratio", numMIR = 1, cutoff = ifelse(method == "ratio", 1/3, 0.05))
ceRNA.basic(miRtar, targetce = NULL, method = "ratio", numMIR = 1, cutoff = ifelse(method == "ratio", 1/3, 0.05))
miRtar |
A data frame representing the relationship between miRNA and target. The data frame contains the name of the miRNA and target regulatory relationship. |
targetce |
a character string (vector) specifying candidate target name to analyse (default (targetce = NULL)). |
method |
a character string (default "ratio") indicating which statistical method to choose to calculate the ceRNA interaction relationship. One of "ratio" (default), or "hypergeometric", can be abbreviated. |
numMIR |
a numeric vector that specify the minimum number of 2 gene-shared miRNAs |
cutoff |
a numeric vector of the method("ratio","hypergeometric") specifying threshold between ceRNA interactions .(default (1/3)). |
Note:All the arguments without default value must be assigned.
A list of identified miRNA sponge interactions containing following components:
cesig
predicted significant triplets and related information,a 5 columns dataframe as following:
targetce
represented target names,respectively.
anotherce
names of modulators that another target(modulators) constitutes a ceRNA interaction relation.
miRNAs
names of miRNA shared by two targets.
miRNAs_num
number of miRNAs shared by two targets.
ratio/pvalue
The ratio/pvalue(optional for method("ratio","hypergeometric")) of the identified ceRNA interaction relation.
cenotsig
predicted not significant triplets and related information,a 5 columns dataframe as following:
targetce
same as targetce
in cesig
anotherce
same as anotherce
in cesig
miRNAs
same as miRNAs
in cesig
miRNAs_num
same as miRNAs_num
in cesig
pvalue
same as pvalue
in cesig
Xu J , Feng L , Han Z , et al. Extensive ceRNA–ceRNA interaction networks mediated by miRNAs regulate development in multiple rhesus tissues[J]. Nucleic Acids Research, 2016:gkw587.
##identifying miRNA sponge interactions ##Here we take six candidate targets(modulators) for example ceRNA.basic(miRtar=dataset[["miRtar"]],targetce=NULL,method="ratio",numMIR = 1,cutoff =1/3)
##identifying miRNA sponge interactions ##Here we take six candidate targets(modulators) for example ceRNA.basic(miRtar=dataset[["miRtar"]],targetce=NULL,method="ratio",numMIR = 1,cutoff =1/3)
identifying miRNA sponge interactions using ceRNA.cernia.In this function We implement cernia score to identify miRNA sponge interactions. The fraction of common miRNAs; 1. The fraction of common miRNAs; 2. The density of the MREs for all shared miRNAs; 3. The distribution of MREs of the putative ceRNAs; 4. The relation between the overall number of MREs for a putative ceRNA, compared with the number of miRNAs that yield these MREs; 5. The density of the hybridization energies related to MREs for all shared miRNAs; 6. The DT-Hybrid recommendation scores; 7. The pairwise Pearson correlation between putative ceRNA expressions.
ceRNA.cernia(miRtar, targetce = NULL, geneexp, miRexp, mres, numMIR = 3, cor_cutoff = 0, s_cutoff = 0.5)
ceRNA.cernia(miRtar, targetce = NULL, geneexp, miRexp, mres, numMIR = 3, cor_cutoff = 0, s_cutoff = 0.5)
miRtar |
A data frame representing the relationship between miRNA and target. The data frame contains the name of the miRNA and target regulatory relationship. |
targetce |
a character string (vector) specifying candidate target name to analyse (default (targetce = NULL)). |
geneexp |
An input target expression data frame, the columns are genes and the rows are samples.The expression value may be gene expression ,non-coding RNA expression or expression values of circRNAs and so on. |
miRexp |
An input miRNA expression data frame, the columns are miRNA and the rows are samples. |
mres |
miRNA response elements (mres) data frame, each row contains five elements: mirna, target, energy, gap_l, gap_r. |
numMIR |
a numeric vector that specify the minimum number of 2 gene-shared miRNAs. |
cor_cutoff |
a numeric vector of the form method specifying threshold between ceRNA interactions correlation,default(cor_cutoff=0). |
s_cutoff |
the threshold of seven comprehensive scores. |
Note:All the arguments without default value must be assigned.The miRNA in the file of the target-miRNA regulatory relationship pair should also be present in the expression profile data file. Internal functions (parMM, graphWeights, recommendation, dtHybrid) of cernia method are from the website: https://github.com/dsardina/cernia Copyright 2016 Rosalba Giugno Licensed under the Apache License, Version 2.0 (the 'License')
A list of identified miRNA sponge interactions containing following components:
targetce
represented target names,respectively.
anotherce
names of modulators that another target(modulators) constitutes a ceRNA interaction relation.
miRNAs_num
names of miRNA in the triplet.
Score 1
The fraction of common miRNAs;
Score 2
The density of the MREs for all shared miRNAs;
Score 3
The distribution of MREs of the putative ceRNAs;
Score 4
The relation between the overall number of MREs for a putative ceRNA,
compared with the number of miRNAs that yield these MREs;
Score 5
The density of the hybridization energies related to MREs for
all shared miRNAs;
Score 6
The DT-Hybrid recommendation scores;
Score 7
The pairwise Pearson correlation between putative ceRNA
expressions.
Normalized score
combine the sum of 7 scores and standardize by z-score.
Sardina D S , Alaimo S , Ferro A , et al. A novel computational method for inferring competing endogenous interactions[J]. Briefings in Bioinformatics, 2016:bbw084.
##identifying miRNA sponge interactions. ##Here we take six candidate targets(modulators) and corresponding expression ##data for example,Specify target(PTEN) to predict ceRNA interaction. ceRNA.cernia(miRtar=dataset[["miRtar"]], targetce = "PTEN", geneexp=dataset[["geneexp"]], numMIR = 1, miRexp=dataset[["miRexp"]], mres=dataset[["mres"]], cor_cutoff = 0.2, s_cutoff = 0.5)
##identifying miRNA sponge interactions. ##Here we take six candidate targets(modulators) and corresponding expression ##data for example,Specify target(PTEN) to predict ceRNA interaction. ceRNA.cernia(miRtar=dataset[["miRtar"]], targetce = "PTEN", geneexp=dataset[["geneexp"]], numMIR = 1, miRexp=dataset[["miRexp"]], mres=dataset[["mres"]], cor_cutoff = 0.2, s_cutoff = 0.5)
identifying miRNA sponge interactions using ceRNA.cmi.In this function We implement CMI methods to identify miRNA sponge interactions.
ceRNA.cmi(miRtar, targetce = NULL, geneexp, miRexp, numMIR = 1, num_perm = 100, cutoff = 0.05)
ceRNA.cmi(miRtar, targetce = NULL, geneexp, miRexp, numMIR = 1, num_perm = 100, cutoff = 0.05)
miRtar |
A data frame representing the relationship between miRNA and target. The data frame contains the name of the miRNA and target regulatory relationship. |
targetce |
a character string (vector) specifying candidate target name to analyse (default (targetce = NULL)). |
geneexp |
An input target expression data frame, the columns are genes and the rows are samples.The expression value may be gene expression ,non-coding RNA expression or expression values of circRNAs and so on. |
miRexp |
An input miRNA expression data frame, the columns are miRNA and the rows are samples. |
numMIR |
a numeric vector that specify the minimum number of 2 gene-shared miRNAs |
num_perm |
The number of random default(num_perm=100). |
cutoff |
a numeric vector of the form method specifying threshold between ceRNA interactions default(cutoff=0.05). |
Note:All the arguments without default value must be assigned.The miRNA in the file of the target-miRNA regulatory relationship pair should also be present in the expression profile data file.
A list of identified miRNA sponge interactions containing following components:
ceRNA_cmi
predicted triplets and related information,a 5 columns dataframe as following:
targetce
represented target names,respectively.
miRNA
names of miRNA in the triplet.
anotherce
names of modulators that another target(modulators) constitutes a ceRNA interaction relation.
cmi
Conditional mutual information(CMI) of triplets calculated using expression values.
pvalue
The p value of the identified ceRNA interaction relation by CMI.
ceRNA_comP
Number of miRNAs interacting with each target in the input file.
targetce
represented target names,respectively.
anotherce
names of predicted another target(modulators) constitutes a ceRNA interaction relation.
commonMiR
names of miRNA shared by predicted ceRNA.
xsq
Chi-square value of p-value pvalue
in conditional mutual information .
comP
The p value calculated by incomplete gamma function .
Sumazin P , Yang X , Chiu H S , et al. An Extensive MicroRNA-Mediated Network of RNA-RNA Interactions Regulates Established Oncogenic Pathways in Glioblastoma[J]. Cell, 2011, 147(2):0-381.
##identifying miRNA sponge interactions. ##Here we take six candidate targets(modulators) and corresponding expression ##data for example,Specify target(PTEN) to predict ceRNA interaction. ceRNA.cmi(miRtar=dataset[["miRtar"]], targetce = "PTEN", geneexp=dataset[["geneexp"]], numMIR = 1, miRexp=dataset[["miRexp"]],num_perm=50)
##identifying miRNA sponge interactions. ##Here we take six candidate targets(modulators) and corresponding expression ##data for example,Specify target(PTEN) to predict ceRNA interaction. ceRNA.cmi(miRtar=dataset[["miRtar"]], targetce = "PTEN", geneexp=dataset[["geneexp"]], numMIR = 1, miRexp=dataset[["miRexp"]],num_perm=50)
identifying miRNA sponge interactions using ceRNA.cor.In this function We implement several popular linear methods ( HyperC, SC) to identify miRNA sponge. interactions.
ceRNA.cor(miRtar, targetce = NULL, geneexp, miRexp, numMIR = 1, method = "pearson", numrandom = 100)
ceRNA.cor(miRtar, targetce = NULL, geneexp, miRexp, numMIR = 1, method = "pearson", numrandom = 100)
miRtar |
A data frame representing the relationship between miRNA and target. The data frame contains the name of the miRNA and target regulatory relationship.Required option for method "pearson" and "partial correlation". |
targetce |
a character string (vector) specifying candidate target name to analyse (default (targetce = NULL)). |
geneexp |
An input target expression data frame, the columns are genes and the rows are samples.The expression value may be gene expression ,non-coding RNA expression or expression values of circRNAs and so on. Required option for method pearson" and "partial correlation". |
miRexp |
An input miRNA expression data frame, the columns are miRNA and the rows are samples. Required option for method pearson" and "partial correlation". |
numMIR |
a numeric vector that specify the minimum number of 2 gene-shared miRNAs. |
method |
a character string (default "pearson") indicating which statistical method to choose to calculate the ceRNA interaction relationship. One of "pearson" (default), or "partial correlation", can be abbreviated. |
numrandom |
The number of random. Required option for method "partial correlation",default (numrandom = 100). |
Note:All the arguments without default value must be assigned.
A list of identified miRNA sponge interactions containing following components:
ceRNA
predicted triplets and related information,a 5 columns dataframe as following:
targetce
represented target names,respectively.
anotherce
names of modulators that another target(modulators) constitutes a ceRNA interaction relation.
miRNAs
names of miRNA shared by two targets.
miRNAs_num
number of miRNAs shared by two targets.
correlation
The correlation of the identified ceRNA interaction relation.
pvalue
The p value of the identified ceRNA interaction relation.
miR_l
Number of miRNAs interacting with each target in the input file.
Paci P , Colombo T , Farina L . Computational analysis identifies a sponge interaction network between long non-coding RNAs and messenger RNAs in human breast cancer[J]. BMC Systems Biology, 2014, 8(1):83. Zhang Y , Xu Y , Feng L , et al. Comprehensive characterization of lncRNA-mRNA related ceRNA network across 12 major cancers[J]. Oncotarget, 2014, 7(39):64148-64167.
##identifying miRNA sponge interactions ##Here we take the regulatory relationship between six genes and 71 miRNAs and corresponding ##expression profilesas an example. ceRNA.cor(miRtar=dataset[["miRtar"]], targetce = NULL, geneexp=dataset[["geneexp"]], miRexp=dataset[["miRexp"]], numMIR = 1,method = "pearson", numrandom = 100)
##identifying miRNA sponge interactions ##Here we take the regulatory relationship between six genes and 71 miRNAs and corresponding ##expression profilesas an example. ceRNA.cor(miRtar=dataset[["miRtar"]], targetce = NULL, geneexp=dataset[["geneexp"]], miRexp=dataset[["miRexp"]], numMIR = 1,method = "pearson", numrandom = 100)
Downstream analysis function of ceRNA,This function can be used to identify biological functions of interest, and users can enrich the functions of interest by ceRNA.enrich.
ceRNA.enrich(data, GOterms, background, threshold = 2, correction = "BH")
ceRNA.enrich(data, GOterms, background, threshold = 2, correction = "BH")
data |
a dataframe that the ceRNA relationship is the data, and the prediction result data obtained according to the ceRNA prediction algorithm.Such as ceRNA.cmi prediction result file. |
GOterms |
a list that goterm of interest and all gene sets in the term. |
background |
a vector containing a gene set in which GOterm annotated genes must be.Its id style must be consistent with the id format in GOterms. |
threshold |
a numeric (default 2) representing min number of intersection between a modulator's targets and a GOterms genes. |
correction |
correction method (default "BH") in one of |
Note:All the arguments without default value must be assigned.
A list of identified miRNA sponge interactions containing following components:
target
represented target names,respectively.
GOterm
the GOterm name.
target_num
names of represented target in the triplet.
GOtermnum
the gene number of a GOterm;
term_tar
the number of intersected factor between a GOterm genes and a modulator targets;
P_value
the p value of the significance enrichment;
fdr
corrected P_value by the assigned method;
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as ceRNA.enrich(data=dataset[["Pre.ceRNA"]],GOterms=dataset[["GOterms"]], background=dataset[["background"]],threshold=1,correction="BH")
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as ceRNA.enrich(data=dataset[["Pre.ceRNA"]],GOterms=dataset[["GOterms"]], background=dataset[["background"]],threshold=1,correction="BH")
A downstream analysis function of ceRNAseekvisualize and analyze the identified ceRNA-ceRNA network, the network can be defined as weighted or un-weighted network.and the basic topological features (such as degree, closeness, betweenness and centrality) of the ceRNAs can be output.
ceRNA.Net(data, net_direct = TRUE, vertex_size = 20, v.label = TRUE, node_shape = "circle",n_color = "orange",E_weight = TRUE, ity = 1, label_cex = 2, label_color = "black",edge_color = "gray", n_frame.color = "gray")
ceRNA.Net(data, net_direct = TRUE, vertex_size = 20, v.label = TRUE, node_shape = "circle",n_color = "orange",E_weight = TRUE, ity = 1, label_cex = 2, label_color = "black",edge_color = "gray", n_frame.color = "gray")
data |
A matrix of ceRNA interaction pairs identified by statistical identification methods selected by the user . |
net_direct |
A logical variable specifies a directed or undirected network.default (net_direct = TRUE). |
vertex_size |
a numeric vector to adjust the node size,default (vertex_size = 20). |
v.label |
Whether to display the label of the node.default (v.label = TRUE). |
node_shape |
A character vector is used to adjust the shape of the node,The selectable shape parameters are "circle","square","csquare","rectangle","crectangle","vrectangle","pie", "sphere","none".default (node_shape = "circle"). |
n_color |
The character vector used to define the fill color of the node. |
E_weight |
A logical vector represents whether the network is a weighted network,default (E_weight = TRUE). |
ity |
A numeric vector defines whether the edge is a solid line or a dotted line,the possible values of the vector are c(1, 2), default (ity = 1). |
label_cex |
Specify the size of the node label font. |
label_color |
Specify the label color of the node. |
edge_color |
Specify the color of the network side. |
n_frame.color |
The character vector used to define the border color of the node. |
This function calls the igraph package. For specific parameter settings, please refer to the igrap help documentation. Note:All the arguments without default value must be assigned.
The output includes two parts: the network diagram of ceRNA interaction and the topology attribute information of the network.
Network topology attributes include 5 types of information:
degree
degree refers to the number of edges in the network directly connected to the node
closeness
An indicator that describes the average distance of a node to all other nodes in the network.
betweenness
The proportion of this node that appears in the shortest path between other nodes.
cluster coefficient
Representing the dense connection nature between some nodes in the network
Eigenvector centrality
Representing the characteristic vector centrality of the network.
##Display ceRNA interactions in a network format and output network topology attributes. ##The input file can be a list of [ceRNA] of the ceRNA.Lin result file, a list of [cesig] ##for the result file identified by ceRNA.basic, or a list of [ceRNA_comP] in the result ##file identified by the ceRNA.cmi function. ##Here we apply the ceRNA list in the example file for CMI identification to ##display the network and analyze the network topology properties. ceRNA.Net(as.matrix(dataset[["Pre.ceRNA"]]),net_direct=TRUE,vertex_size=20,v.label = TRUE, node_shape="circle",n_color = "orange",E_weight=TRUE,ity=1,label_cex=2, label_color="black",edge_color="gray",n_frame.color = "gray")
##Display ceRNA interactions in a network format and output network topology attributes. ##The input file can be a list of [ceRNA] of the ceRNA.Lin result file, a list of [cesig] ##for the result file identified by ceRNA.basic, or a list of [ceRNA_comP] in the result ##file identified by the ceRNA.cmi function. ##Here we apply the ceRNA list in the example file for CMI identification to ##display the network and analyze the network topology properties. ceRNA.Net(as.matrix(dataset[["Pre.ceRNA"]]),net_direct=TRUE,vertex_size=20,v.label = TRUE, node_shape="circle",n_color = "orange",E_weight=TRUE,ity=1,label_cex=2, label_color="black",edge_color="gray",n_frame.color = "gray")
It is used to predict the survival of ternary relationship pairs and to support the survival prognosis of training sets and test sets.
ceRNA.surv(ceRNA, exp.sur, train = NULL, test = NULL, index)
ceRNA.surv(ceRNA, exp.sur, train = NULL, test = NULL, index)
ceRNA |
a dataframe that the ceRNA relationship is the data, and the prediction result data obtained according to the ceRNA prediction algorithm.Such as ceRNA.cmi prediction result file. |
exp.sur |
dataframe specifying expression and survival information.Its rownames are sample names.Its colnames are names in triplets and survival information (see example data in details). |
train |
a character string vector specifying train sample names. |
test |
a character string vector specifying test sample names. |
index |
a numeric vector (default 1) reprsenting rowindex of triplets analyed. |
A list of identified miRNA sponge interactions containing following components:
targetce
represented target names,respectively.
anotherce
names of modulators that another target(modulators) constitutes a ceRNA interaction relation.
coef_targetce
the coxph coefficient of targetce;
p_targetce
the coxph significance of targetce;
coef_anotherce
the coxph coefficient of anotherce;
p_anotherce
the coxph significance of anotherce;
N_low
The genes were ranked according to expression, and the number of samples
expressed in the bottom 25
N_high
The genes were ranked according to expression, and the number of samples
expressed in the top 25
HR
Risk score of ceRNA ternary pair.
p
Survival significance of ceRNA ternary pairs.
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as ceRNA.surv(ceRNA=dataset[["Pre.ceRNA"]],exp.sur=dataset[["exp.sur"]],train=NULL,test=NULL,index=5)
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as ceRNA.surv(ceRNA=dataset[["Pre.ceRNA"]],exp.sur=dataset[["exp.sur"]],train=NULL,test=NULL,index=5)
This object contains data for examples.
data(dataset)
data(dataset)
A list with 4 variables:
MiRNA–target regulation,assembled genome-wide miRNA-gene regulation by TargetScan,It includes 44 rows standing for miRNA-gene regulation pairs in glioblastoma and 2 columns (the first is dataframe header).
a dataframe representing target regulations expression ,all of which are gene expression value,in glioblastoma.It includes 6 rows standing for 6 gene targets and 541 columns standing for 541 samples.Its rownames are gene symbols.
a dataframe representing miRNA expression profile in glioblastoma.It includes 57 rows standing for 57 miRNA and 541 columns standing for 541 samples.Its rownames are miRNA symbols.
a data frame representing the ceRNA identified in the pre-experiment using the method provided by the software to draw a network map.
a data frame representing the MRE.
a list representing that goterm of interest and all gene sets in the term.
a character representing a gene set in which GOterm annotated genes must be.
a data frame representing specifying expression and survival information.Its rownames are sample names.Its colnames are names in triplets and survival information
a data frame representing specifying train sample names
a data frame representing specifying test sample names.
All epression data is from a study about glioblastoma.The miRNA-gene regulation information is from TargetScan.
Network, Atlas T C G. et al.
Network, Atlas T C G . Corrigendum: Comprehensive genomic characterization defines human glioblastoma genes and core pathways[J]. Nature, 2013, 494(7438):506-506.
data(dataset) ## maybe str(dataset) ; ...
data(dataset) ## maybe str(dataset) ; ...
It is used to draw the survival curve of ternary relationship pairs and to support the survival prognosis of training sets and test sets.
surv.plot(ceRNA, exp.sur, train = NULL, test = NULL, index)
surv.plot(ceRNA, exp.sur, train = NULL, test = NULL, index)
ceRNA |
a dataframe that the ceRNA relationship is the data, and the prediction result data obtained according to the ceRNA prediction algorithm.Such as ceRNA.cmi prediction result file. |
exp.sur |
dataframe specifying expression and survival information.Its rownames are sample names.Its colnames are names in triplets and survival information (see example data in details). |
train |
a character string vector specifying train sample names. |
test |
a character string vector specifying test sample names. |
index |
a numeric vector (default 1) reprsenting rowindex of triplets analyed. |
Survival curve of user-selected ceRNA ternary pairs.
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as surv.plot(ceRNA=dataset[["Pre.ceRNA"]],exp.sur=dataset[["exp.sur"]],train=NULL,test=NULL,index=5)
##---- Should be DIRECTLY executable !! ---- ##-- ==> Define data, use random, ##-- or do help(data=index) for the standard data sets. ## The function is currently defined as surv.plot(ceRNA=dataset[["Pre.ceRNA"]],exp.sur=dataset[["exp.sur"]],train=NULL,test=NULL,index=5)