Title: | Gene Regulatory Network Inference Using Time Series |
---|---|
Description: | The package offers four network inference statistical models using Dynamic Bayesian Networks and Gibbs Variable Selection: a linear interaction model, two linear interaction models with added experimental noise (Gaussian and Student distributed) for the case where replicates are available and a non-linear interaction model. |
Authors: | Edward Morrissey |
Maintainer: | Edward Morrissey <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.59.0 |
Built: | 2024-12-21 05:57:20 UTC |
Source: | https://github.com/bioc/GRENITS |
Analyse output from network inference functions. Basic convergence and analysis plots.
analyse.output(output.folder, timeSeries = NULL)
analyse.output(output.folder, timeSeries = NULL)
output.folder |
Name of folder (including path) where chains are kept |
timeSeries |
Only used by NonLinearModel analysis. Data matrix containing gene expression time series. Where genes will be placed in rows and time points in columns. |
Read first two chains run and plot some basic convergence plots (ConvergencePlots.pdf), analysis plots (AnalysisPlots.pdf), as well as inferred network probabilities in two formats (NetworkProbability_List.txt and NetworkProbability_Matrix.txt).
The output of the analysis will be four files (five if nonLinearNet). The contents of the two plot files change depending on the network inference function used.
ConvergencePlots.pdf |
Basic convergence plots. The posterior means of each variable are compared. |
AnalysisPlots.pdf |
Heatmap plot of network link probabilities as well as marginal network uncertainty plot. A plot of the number of links predicted by the model for a given probability threshold. For ReplicatesNet_student, the posterior distribution of the degrees of freedom are also plotted. For NonLinearNet, the posterior of the smoothness parameter is plotted. |
NetworkProbability_List.txt |
Posterior probabilities for each network connection in list format, including posterior interaction strength for linear models. Can be imported with network analysis software such as cytoscape. |
NetworkProbability_Matrix.txt |
Posterior probabilities for each network connection in matrix format. |
ProbNumParents.txt |
Posterior probabilities for number of regulators for each gene. |
InferredFunctionPlots.pdf |
(Only for nonLinearNet) Posterior distribution of predicted functions. Data values are plotted as circles. |
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2010. On reverse engineering of gene interaction networks using time course data with repeated measurements. Bioinformatics 2010; doi: 10.1093/bioinformatics/btq421
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2011 Inferring the time-invariant topology of a nonlinear sparse gene regulatory network using fully Bayesian spline autoregression Biostatistics 2011; doi: 10.1093/biostatistics/kxr009
NonLinearNet
, LinearNet
,
ReplicatesNet_student
, ReplicatesNet_gauss
.
# Load A. thaliana circadian clock ODE generated data data(Athaliana_ODE) # Folder where raw runs will be kept and analysed output.folder <- paste(tempdir(), "/Example_LinearNet",sep="") # Run network inference, place raw results in output.folder LinearNet(output.folder, Athaliana_ODE) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder)
# Load A. thaliana circadian clock ODE generated data data(Athaliana_ODE) # Folder where raw runs will be kept and analysed output.folder <- paste(tempdir(), "/Example_LinearNet",sep="") # Run network inference, place raw results in output.folder LinearNet(output.folder, Athaliana_ODE) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder)
This data set was generated using a circadian clock regulatory network ODE model. The data was simulated using COPASI, subsampled to produce hourly data and after this the log was calculated.
data(Athaliana_ODE)
data(Athaliana_ODE)
A matrix containing genes in rows and (1h) time points in columns.
Hoops, S., Sahle, S., Gauges, R., Lee, C., Pahle, J., Simus, N., Singhal, M., Xu, L., Mendes, P. and Kummer, U. (2006) COPASI a COmplex PAthway SImulator. Bioinformatics, 22, 3067-3074.
Locke, J.C.W., Kozma-Bognar, L., Gould, P.D., Feher, B., Kevei, E., Nagy, F., Turner, M.S., Hall, A. and Millar, A.J. (2006) Experimental validation of a predicted feedback loop in the multi-oscillator clock of Arabidopsis thaliana. Molecular Systems Biology
This data set was generated using a circadian clock regulatory network ODE model after which noise was added. The data was simulated using COPASI, subsampled to produce hourly data and after this the log was caluclated. Further to this, four replicates were produced by adding student noise.
data(Athaliana_ODE_4NoiseReps)
data(Athaliana_ODE_4NoiseReps)
A matrix containing genes in rows and (1h) time points in columns. The four replicates are appended (columns)
Hoops, S., Sahle, S., Gauges, R., Lee, C., Pahle, J., Simus, N., Singhal, M., Xu, L., Mendes, P. and Kummer, U. (2006) COPASI a COmplex PAthway SImulator. Bioinformatics, 22, 3067-3074.
Locke, J.C.W., Kozma-Bognar, L., Gould, P.D., Feher, B., Kevei, E., Nagy, F., Turner, M.S., Hall, A. and Millar, A.J. (2006) Experimental validation of a predicted feedback loop in the multi-oscillator clock of Arabidopsis thaliana. Molecular Systems Biology
Run Bayesian inference of linear interaction network. The function generates MCMC chains that can later be analysed.
LinearNet( resultsFolder, timeSeries, ParamVec = NULL, chains = 2, user.seeds = NULL, Regulators = NULL, fixMe = NULL)
LinearNet( resultsFolder, timeSeries, ParamVec = NULL, chains = 2, user.seeds = NULL, Regulators = NULL, fixMe = NULL)
resultsFolder |
Name of output folder. The folder will be created and the output of the run will be placed there. |
timeSeries |
Data matrix containing gene expression time series. Where genes will be placed in rows and time points in columns. Gene names may be included as row names. |
ParamVec |
A parameter vector created using "mcmc.defaultParams_Linear". If none is given, default parameters will be used. The vector contains parameters associated to the priors as well as MCMC run length. (See mcmc.defaultParams_Linear) |
chains |
Number of MCMC chains to run. |
user.seeds |
An optional vector with seeds to use for MCMC chains. |
Regulators |
An optional vector with the indices of which genes are regulators. If provided, all non-regulator genes will not be allowed to regulate. |
fixMe |
An optional matrix of size genes x genes, where columns represent regulators and rows regulated genes. The matrix informs the model of network connections known to be present/absent. For each position use either 0 (no regulation, fix off), 1 (known regulatory interaction, fix on) or NaN (no information, do not fix). |
For each chain run, a folder (chain1, chain2, ...) will be created and the output of the MCMC run will be placed there. The files will be B_mcmc (the coeffcients of the linear regression), Gamma_mcmc (the indicator variables of Gibbs variable selection), Lambda_mcmc (the precision of each regression), Mu_mcmc (the intercept of each regression) and Rho_mcmc (the network connectivity parameter).
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2010. On reverse engineering of gene interaction networks using time course data with repeated measurements. Bioinformatics 2010; doi: 10.1093/bioinformatics/btq421
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2011 Inferring the time-invariant topology of a nonlinear sparse gene regulatory network using fully Bayesian spline autoregression Biostatistics 2011; doi: 10.1093/biostatistics/kxr009
mcmc.defaultParams_Linear, analyse.output
.
# Load A. thaliana circadian clock ODE generated data data(Athaliana_ODE) # Folder where raw runs will be kept and analysed output.folder <- paste(tempdir(), "/Example_LinearNet",sep="") # Run network inference, place raw results in output.folder LinearNet(output.folder, Athaliana_ODE) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder)
# Load A. thaliana circadian clock ODE generated data data(Athaliana_ODE) # Folder where raw runs will be kept and analysed output.folder <- paste(tempdir(), "/Example_LinearNet",sep="") # Run network inference, place raw results in output.folder LinearNet(output.folder, Athaliana_ODE) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder)
Create parameter vector with default parameters for ReplicatesNet_gauss function
mcmc.defaultParams_gauss()
mcmc.defaultParams_gauss()
Use this function to generate a template parameter vector to use non-default parameters for the ReplicatesNet_gauss model.
Returns a single vector with the following elements (in this order):
(1) samples |
Number of MCMC iterations to run |
(2) burn.in |
Number of initial iterations to discard as burn in |
(3) thin |
Subsampling frequency |
(4) c |
Shape parameter 1 for Beta(c,d) prior on rho (connectivity parameter) |
(5) d |
Shape parameter 2 for Beta(c,d) prior on rho (connectivity parameter) |
(6) sigma.s |
Standard deviation parameter for N(0,sigma.s) prior on B (Regression coefficients) |
(7) a |
Shape parameter for Gamma(a,b) prior on lambda (Regression precision) |
(8) b |
Rate parameter for Gamma(a,b) prior on lambda (Regression precision) |
(9) a_exp |
Shape parameter for Gamma(a_exp,b_exp) prior on tau (Replicates precision) |
(10) b_exp |
Rate parameter for Gamma(a_exp,b_exp) prior on tau (Replicates precision) |
(11) sigma.mu |
Standard deviation parameter for N(0,sigma.mu) prior on mu (Regression intercept) |
(12) fix.y.iter |
Number of iterations for which sampled data Y is fixed |
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2010. On reverse engineering of gene interaction networks using time course data with repeated measurements. Bioinformatics 2010; doi: 10.1093/bioinformatics/btq421
plotPriors
, ReplicatesNet_gauss
.
# Get default parameters linearNet_Gauss.params <- mcmc.defaultParams_gauss() # Change run length linearNet_Gauss.params[1] <- 200000 # Change prior regression precision linearNet_Gauss.params[7] <- 0.001 linearNet_Gauss.params[8] <- 0.001 # Plot to visualise changes plotPriors(linearNet_Gauss.params) ## Use to run ReplicatesNet_gauss ...
# Get default parameters linearNet_Gauss.params <- mcmc.defaultParams_gauss() # Change run length linearNet_Gauss.params[1] <- 200000 # Change prior regression precision linearNet_Gauss.params[7] <- 0.001 linearNet_Gauss.params[8] <- 0.001 # Plot to visualise changes plotPriors(linearNet_Gauss.params) ## Use to run ReplicatesNet_gauss ...
Create parameter vector with default parameters for LinearNet function
mcmc.defaultParams_Linear()
mcmc.defaultParams_Linear()
Use this function to generate a template parameter vector to use non-default parameters for the LinearNet model.
Returns a single vector with the following elements (in this order):
(1) samples |
Number of MCMC iterations to run. |
(2) burn.in |
Number of initial iterations to discard as burn in. |
(3) thin |
Subsampling frequency |
(4) c |
Shape parameter 1 for Beta(c,d) prior on rho (connectivity parameter) |
(5) d |
Shape parameter 2 for Beta(c,d) prior on rho (connectivity parameter) |
(6) sigma.s |
Standard deviation parameter for N(0,sigma.s) prior on B (Regression coefficients) |
(7) a |
Shape parameter for Gamma(a,b) prior on lambda (Regression precision) |
(8) b |
Rate parameter for Gamma(a,b) prior on lambda (Regression precision) |
(9) sigma.mu |
Standard deviation parameter for N(0,sigma.mu) prior on mu (Regression intercept) |
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2010. On reverse engineering of gene interaction networks using time course data with repeated measurements. Bioinformatics 2010; doi: 10.1093/bioinformatics/btq421
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2011 Inferring the time-invariant topology of a nonlinear sparse gene regulatory network using fully Bayesian spline autoregression Biostatistics 2011; doi: 10.1093/biostatistics/kxr009
# Get default parameters linearNet.params <- mcmc.defaultParams_Linear() # Change run length linearNet.params[1] <- 150000 # Change prior regression precision linearNet.params[7] <- 0.001 linearNet.params[8] <- 0.001 # Plot to check changes plotPriors(linearNet.params) ## Use to run LinearNet ...
# Get default parameters linearNet.params <- mcmc.defaultParams_Linear() # Change run length linearNet.params[1] <- 150000 # Change prior regression precision linearNet.params[7] <- 0.001 linearNet.params[8] <- 0.001 # Plot to check changes plotPriors(linearNet.params) ## Use to run LinearNet ...
Create parameter vector with default parameters for NonLinearNet function
mcmc.defaultParams_nonLinear()
mcmc.defaultParams_nonLinear()
Use this function to generate a template parameter vector to use non-default parameters for the NonLinearNet model.
Returns a single vector with the following elements (in this order):
(1) samples |
Number of MCMC iterations to run. |
(2) burn.in |
Number of initial iterations to discard as burn in. |
(3) thin |
Subsampling frequency |
(4) c |
Shape parameter 1 for Beta(c,d) prior on rho (connectivity parameter) |
(5) d |
Shape parameter 2 for Beta(c,d) prior on rho (connectivity parameter) |
(6) trunc |
Truncation parameter for InvertedPareto prior on tau (smoothness parameter) |
(7) tau0 |
Precision parameter for N(0, tau0^(-0.5)) prior on B (first two coefficients) |
(8) M |
Numer of knots used for each spline function |
(9) a |
Shape parameter for Gamma(a,b) prior on lambda (Regression precision) |
(10) b |
Rate parameter for Gamma(a,b) prior on lambda (Regression precision) |
(11) sigma.mu |
Standard deviation parameter for N(0,sigma.mu) prior on mu (Regression intercept) |
(12) a_pareto |
Pareto parameter for InvertedPareto prior on tau (smoothness parameter) |
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2011 Inferring the time-invariant topology of a nonlinear sparse gene regulatory network using fully Bayesian spline autoregression Biostatistics 2011; doi: 10.1093/biostatistics/kxr009
# Get default parameters nonLinearNet.params <- mcmc.defaultParams_nonLinear() # Change run length nonLinearNet.params[1] <- 150000 # Change prior on smoothness parameter nonLinearNet.params[6] <- 30000 # Change truncation nonLinearNet.params[12] <- 3 # Concentrate more mass close to linear region # Plot to check changes plotPriors(nonLinearNet.params) ## Use to run LinearNet ...
# Get default parameters nonLinearNet.params <- mcmc.defaultParams_nonLinear() # Change run length nonLinearNet.params[1] <- 150000 # Change prior on smoothness parameter nonLinearNet.params[6] <- 30000 # Change truncation nonLinearNet.params[12] <- 3 # Concentrate more mass close to linear region # Plot to check changes plotPriors(nonLinearNet.params) ## Use to run LinearNet ...
Create parameter vector with default parameters for ReplicatesNet_student function
mcmc.defaultParams_student()
mcmc.defaultParams_student()
Use this function to generate a template parameter vector to use non-default parameters for the ReplicatesNet_student model.
Returns a single vector with the following elements (in this order):
(1) samples |
Number of MCMC iterations to run |
(2) burn.in |
Number of initial iterations to discard as burn in |
(3) thin |
Subsampling frequency |
(4) c |
Shape parameter 1 for Beta(c,d) prior on rho (connectivity parameter) |
(5) d |
Shape parameter 2 for Beta(c,d) prior on rho (connectivity parameter) |
(6) sigma.s |
Standard deviation parameter for N(0,sigma.s) prior on B (Regression coefficients) |
(7) a |
Shape parameter for Gamma(a,b) prior on lambda (Regression precision) |
(8) b |
Rate parameter for Gamma(a,b) prior on lambda (Regression precision) |
(9) a_exp |
Shape parameter for Gamma(a_exp,b_exp) prior on tau (Replicates precision) |
(10) b_exp |
Rate parameter for Gamma(a_exp,b_exp) prior on tau (Replicates precision) |
(11) a_deg |
Shape parameter for Gamma(a_deg,b_deg) prior on nu (Student degrees of freedom) |
(12) b_deg |
Rate parameter for Gamma(a_deg,b_deg) prior on nu (Student degrees of freedom) |
(13) sigma.mu |
Standard deviation parameter for N(0,sigma.mu) prior on mu (Regression intercept) |
(14) fix.y.iter |
Number of iterations for which sampled data Y is fixed |
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2010. On reverse engineering of gene interaction networks using time course data with repeated measurements. Bioinformatics 2010; doi: 10.1093/bioinformatics/btq421
plotPriors
, ReplicatesNet_student
.
# Get default parameters linearNet_Student.params <- mcmc.defaultParams_student() # Change run length linearNet_Student.params[1] <- 200000 # Change prior regression precision linearNet_Student.params[7] <- 0.001 linearNet_Student.params[8] <- 0.001 # Plot to visualise changes plotPriors(linearNet_Student.params) ## Use to run ReplicatesNet_student ...
# Get default parameters linearNet_Student.params <- mcmc.defaultParams_student() # Change run length linearNet_Student.params[1] <- 200000 # Change prior regression precision linearNet_Student.params[7] <- 0.001 linearNet_Student.params[8] <- 0.001 # Plot to visualise changes plotPriors(linearNet_Student.params) ## Use to run ReplicatesNet_student ...
Run Bayesian inference of non-linear interaction network. Non linear interactions are modelled using Penalised Splines. The function generates MCMC chains that can later be analysed.
NonLinearNet( resultsFolder, timeSeries, ParamVec = NULL, chains = 2, user.seeds = NULL, Regulators = NULL, fixMe = NULL)
NonLinearNet( resultsFolder, timeSeries, ParamVec = NULL, chains = 2, user.seeds = NULL, Regulators = NULL, fixMe = NULL)
resultsFolder |
Name of output folder. The folder will be created and the output of the run will be placed there. |
timeSeries |
Data matrix containing gene expression time series. Where genes will be placed in rows and time points in columns. Gene names may be included as row names. |
ParamVec |
A parameter vector created using "mcmc.defaultParams_nonLinear". If none is given, default parameters will be used. The vector contains parameters associated to the priors as well as MCMC run length. (See mcmc.defaultParams_nonLinear) |
chains |
Number of MCMC chains to run. |
user.seeds |
An optional vector with seeds to use for MCMC chains. |
Regulators |
An optional vector with the indices of which genes are regulators. If provided, all non-regulator genes will not be allowed to regulate. |
fixMe |
An optional matrix of size genes x genes, where columns represent regulators and rows regulated genes. The matrix informs the model of network connections known to be present/absent. For each position use either 0 (no regulation, fix off), 1 (known regulatory interaction, fix on) or NaN (no information, do not fix). |
For each chain run, a folder (chain1, chain2, ...) will be created and the output of the MCMC run will be placed there. The files will be Gamma_mcmc (the indicator variables of Gibbs variable selection), Lambda_mcmc (the precision of each regression), Mu_mcmc (the intercept of each regression), Rho_mcmc (the network connectivity parameter), Tau_mcmc (the "smoothness parameter"), all_f (posterior mean of all functions), all_f_sqr (posterior mean of the square of all functions) and Full_F_sqr (posterior mean of the square of the sum of all functions, for each regression). For the files all_f and all_f_sqr functions are placed in column-wise order. The file is filled by placing all interactions for each regression one after another.
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2011 Inferring the time-invariant topology of a nonlinear sparse gene regulatory network using fully Bayesian spline autoregression Biostatistics 2011; doi: 10.1093/biostatistics/kxr009
mcmc.defaultParams_nonLinear, analyse.output
.
# Synthetic data data(Athaliana_ODE) # Reduced data set to 3 genes 20 TP for faster run Athaliana_ODE.reduced <- Athaliana_ODE[c(1,3,5),1:20] # Folder where raw runs will be kept and later analysed output.folder <- paste(tempdir(), "/ExampleNonLinearNet", sep = "") # Run network inference and place raw results in output.folder NonLinearNet(output.folder , Athaliana_ODE.reduced) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder, Athaliana_ODE.reduced)
# Synthetic data data(Athaliana_ODE) # Reduced data set to 3 genes 20 TP for faster run Athaliana_ODE.reduced <- Athaliana_ODE[c(1,3,5),1:20] # Folder where raw runs will be kept and later analysed output.folder <- paste(tempdir(), "/ExampleNonLinearNet", sep = "") # Run network inference and place raw results in output.folder NonLinearNet(output.folder , Athaliana_ODE.reduced) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder, Athaliana_ODE.reduced)
Plot appropriate priors using parameters from vector
plotPriors(parameter.vec)
plotPriors(parameter.vec)
parameter.vec |
MCMC parameter vector of the type generated by e.g. mcmc.defaultParams_Linear |
This function takes the parameter vector that will be used for network inference function and plots the priors associated with the parameters given.
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2010. On reverse engineering of gene interaction networks using time course data with repeated measurements. Bioinformatics 2010; doi: 10.1093/bioinformatics/btq421
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2011 Inferring the time-invariant topology of a nonlinear sparse gene regulatory network using fully Bayesian spline autoregression Biostatistics 2011; doi: 10.1093/biostatistics/kxr009
mcmc.defaultParams_gauss
, mcmc.defaultParams_Linear
,
mcmc.defaultParams_nonLinear
, mcmc.defaultParams_student
.
# Get default parameters nonLinearNet.params <- mcmc.defaultParams_nonLinear() # Change run length nonLinearNet.params[1] <- 150000 # Change prior on smoothness parameter nonLinearNet.params[6] <- 30000 # Change truncation nonLinearNet.params[12] <- 3 # Concentrate more mass close to linear region # Plot to check changes plotPriors(nonLinearNet.params)
# Get default parameters nonLinearNet.params <- mcmc.defaultParams_nonLinear() # Change run length nonLinearNet.params[1] <- 150000 # Change prior on smoothness parameter nonLinearNet.params[6] <- 30000 # Change truncation nonLinearNet.params[12] <- 3 # Concentrate more mass close to linear region # Plot to check changes plotPriors(nonLinearNet.params)
Read MCMC chains for further analysis.
read.chain(output.folder, chainNumber)
read.chain(output.folder, chainNumber)
output.folder |
Name of folder (including path) where chains are kept |
chainNumber |
Which of the chains will be read |
Read chains produced by NonLinearNet, LinearNet, ReplicatesNet_student and ReplicatesNet_gauss for further analysis.
Returns a list of vectors/matrices with the value of the variables at each MCMC iteration.
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2010. On reverse engineering of gene interaction networks using time course data with repeated measurements. Bioinformatics 2010; doi: 10.1093/bioinformatics/btq421
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2011 Inferring the time-invariant topology of a nonlinear sparse gene regulatory network using fully Bayesian spline autoregression Biostatistics 2011; doi: 10.1093/biostatistics/kxr009
NonLinearNet
, LinearNet
,
ReplicatesNet_student
, ReplicatesNet_gauss
.
############################################# ## Run inference using one chain ############################################# # Load A. thaliana circadian clock ODE generated data data(Athaliana_ODE) # Folder where raw runs will be kept and analysed output.folder <- paste(tempdir(), "/Example_LinearNet", sep="") # Run network inference, place raw results in output.folder # Run just one chain for example purpose LinearNet(output.folder, Athaliana_ODE, chains = 1) ########################### ## Read chain ########################### chain1 <- read.chain(output.folder, 1)
############################################# ## Run inference using one chain ############################################# # Load A. thaliana circadian clock ODE generated data data(Athaliana_ODE) # Folder where raw runs will be kept and analysed output.folder <- paste(tempdir(), "/Example_LinearNet", sep="") # Run network inference, place raw results in output.folder # Run just one chain for example purpose LinearNet(output.folder, Athaliana_ODE, chains = 1) ########################### ## Read chain ########################### chain1 <- read.chain(output.folder, 1)
Run Bayesian inference of linear interaction network on data with replicates. The replicates are assumed to follow a Gaussian distribution. The function generates MCMC chains that can later be analysed.
ReplicatesNet_gauss( resultsFolder, timeSeries, numReps, ParamVec = NULL, chains = 2, user.seeds = NULL, Regulators = NULL, fixMe = NULL)
ReplicatesNet_gauss( resultsFolder, timeSeries, numReps, ParamVec = NULL, chains = 2, user.seeds = NULL, Regulators = NULL, fixMe = NULL)
resultsFolder |
Name of output folder. The folder will be created and the output of the run will be placed there. |
timeSeries |
Data matrix containing gene expression time series. Where genes will be placed in rows and time points in columns. Each times series must be placed one after another, so that the first n columns correspond to time series replicate one, the next n columns to time series replicate two, etc. Gene names may be included as row names. |
numReps |
Number of replicate time series included in timeSeries matrix. |
ParamVec |
A parameter vector created using "mcmc.defaultParams_gauss". If none is given, default parameters will be used. The vector contains parameters associated to the priors as well as MCMC run length. (See mcmc.defaultParams_gauss) |
chains |
Number of MCMC chains to run. |
user.seeds |
An optional vector with seeds to use for MCMC chains. |
Regulators |
An optional vector with the indices of which genes are regulators. If provided, all non-regulator genes will not be allowed to regulate. |
fixMe |
An optional matrix of size genes x genes, where columns represent regulators and rows regulated genes. The matrix informs the model of network connections known to be present/absent. For each position use either 0 (no regulation, fix off), 1 (known regulatory interaction, fix on) or NaN (no information, do not fix). |
The order in which the replicates are placed do not affect the output. In other words swapping timepoint 2 replicate 1 and timepoint 2 replicate 2 makes no difference. For the cases where a measurement is not available, an NaN may be used.
For each chain run, a folder (chain1, chain2, ...) will be created and the output of the MCMC run will be placed there. The files will be B_mcmc (the coeffcients of the linear regression), Gamma_mcmc (the indicator variables of Gibbs variable selection), Lambda_mcmc (the precision of each regression), Mu_mcmc (the intercept of each regression), Rho_mcmc (the network connectivity parameter), DataMean_standarised (a times series with the mean of the replicates), Lambda_exp_mcmc (the precision of the replicate noise) and Y_mean (the mean posterior value of the inferred "true mRNA").
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2010. On reverse engineering of gene interaction networks using time course data with repeated measurements. Bioinformatics 2010; doi: 10.1093/bioinformatics/btq421
mcmc.defaultParams_gauss
, analyse.output
.
# Synthetic data data(Athaliana_ODE_4NoiseReps) # Folder where raw runs will be kept and later analysed output.folder <- paste(tempdir(), "/Example_ReplicatesNet_Gauss",sep="") # Run network inference and place raw results in output.folder ReplicatesNet_gauss(output.folder, Athaliana_ODE_4NoiseReps, numReps = 4) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder)
# Synthetic data data(Athaliana_ODE_4NoiseReps) # Folder where raw runs will be kept and later analysed output.folder <- paste(tempdir(), "/Example_ReplicatesNet_Gauss",sep="") # Run network inference and place raw results in output.folder ReplicatesNet_gauss(output.folder, Athaliana_ODE_4NoiseReps, numReps = 4) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder)
Run Bayesian inference of linear interaction network on data with replicates. The replicates are assumed to follow a Student distribution. The function generates MCMC chains that can later be analysed.
ReplicatesNet_student( resultsFolder, timeSeries, numReps, ParamVec = NULL, chains = 2, user.seeds = NULL, Regulators = NULL, fixMe = NULL)
ReplicatesNet_student( resultsFolder, timeSeries, numReps, ParamVec = NULL, chains = 2, user.seeds = NULL, Regulators = NULL, fixMe = NULL)
resultsFolder |
Name of output folder. The folder will be created and the output of the run will be placed there. |
timeSeries |
Data matrix containing gene expression time series. Where genes will be placed in rows and time points in columns. Each times series must be placed one after another, so that the first n columns correspond to time series replicate one, the next n columns to time series replicate two, etc. Gene names may be included as row names. |
numReps |
Number of replicate time series included in timeSeries matrix. |
ParamVec |
A parameter vector created using "mcmc.defaultParams_student". If none is given, default parameters will be used. The vector contains parameters associated to the priors as well as MCMC run length. (See mcmc.defaultParams_student) |
chains |
Number of MCMC chains to run. |
user.seeds |
An optional vector with seeds to use for MCMC chains. |
Regulators |
An optional vector with the indices of which genes are regulators. If provided, all non-regulator genes will not be allowed to regulate. |
fixMe |
An optional matrix of size genes x genes, where columns represent regulators and rows regulated genes. The matrix informs the model of network connections known to be present/absent. For each position use either 0 (no regulation, fix off), 1 (known regulatory interaction, fix on) or NaN (no information, do not fix). |
The order in which the replicates are placed do not affect the output. In other words swapping timepoint 2 replicate 1 and timepoint 2 replicate 2 makes no difference. For the cases where a measurement is not available, an NaN may be used.
For each chain run, a folder (chain1, chain2, ...) will be created and the output of the MCMC run will be placed there. The files will be B_mcmc (the coeffcients of the linear regression), Gamma_mcmc (the indicator variables of Gibbs variable selection), Lambda_mcmc (the precision of each regression), Mu_mcmc (the intercept of each regression), Rho_mcmc (the network connectivity parameter), DataMean_standarised (a times series with the mean of the replicates), Lambda_exp_mcmc (the precision of the replicate noise), Y_mean (the mean posterior value of the inferred "true mRNA"), DegFreedom_mcmc (degrees of freedom of the student distributed replicate noise) and acceptanceRatio (acceptance ratio of the Metropolis-Hastings degrees of freedom update).
Morrissey, E.R., Juarez, M.A., Denby, K.J. and Burroughs, N.J. 2010. On reverse engineering of gene interaction networks using time course data with repeated measurements. Bioinformatics 2010; doi: 10.1093/bioinformatics/btq421
mcmc.defaultParams_student
, analyse.output
.
# Synthetic data data(Athaliana_ODE_4NoiseReps) # Folder where raw runs will be kept and later analysed output.folder <- paste(tempdir(), "/Example_ReplicatesNet_Student",sep="") # Run network inference and place raw results in output.folder ReplicatesNet_student(output.folder, Athaliana_ODE_4NoiseReps, numReps = 4) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder)
# Synthetic data data(Athaliana_ODE_4NoiseReps) # Folder where raw runs will be kept and later analysed output.folder <- paste(tempdir(), "/Example_ReplicatesNet_Student",sep="") # Run network inference and place raw results in output.folder ReplicatesNet_student(output.folder, Athaliana_ODE_4NoiseReps, numReps = 4) # Analyse raw results, place analysis plots and files in output.folder analyse.output(output.folder)