Title: | Dynamic Transcriptome Analysis |
---|---|
Description: | Dynamic Transcriptome Analysis (DTA) can monitor the cellular response to perturbations with higher sensitivity and temporal resolution than standard transcriptomics. The package implements the underlying kinetic modeling approach capable of the precise determination of synthesis- and decay rates from individual microarray or RNAseq measurements. |
Authors: | Bjoern Schwalb, Benedikt Zacher, Sebastian Duemcke, Achim Tresch |
Maintainer: | Bjoern Schwalb <[email protected]> |
License: | Artistic-2.0 |
Version: | 2.53.0 |
Built: | 2024-11-19 03:27:17 UTC |
Source: | https://github.com/bioc/DTA |
The DTA
package implements all methods of the quantitative kinetic modeling approach belonging to DTA (Dynamic Transcriptome Analysis) to estimate mRNA synthesis and decay rates from individual time point measurements.
Package: | DTA |
Type: | Package |
Version: | 2.0.1 |
Date: | 2012-03-22 |
License: | Artistic-2.0 |
LazyLoad: | yes |
Bjoern Schwalb <[email protected]>
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011. M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review. B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics.
## see vignette or supplemental material of the given references.
## see vignette or supplemental material of the given references.
The amount of thymines in the cDNA of each transcript of all Drosophila Melanogaster Ensembl transcript IDs (Flybase transcript number), to assess the uridine-dependent labeling bias and eventually correct for it.
Dm.tnumber
Dm.tnumber
Vector gives the number of thymines in the cDNA (uridine residues in RNA) of each Ensembl transcript ID.
E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
R object contains all relevant *.RData files needed for the DTA.estimate
function. For example, see vignette.
Doelken2008
Doelken2008
R object contains the following *.RData files:
Hs.phenomat
Hs.datamat
Hs.reliable
Hs.enst2ensg
Hs.tnumber
Mm.phenomat
Mm.datamat
Mm.reliable
Mm.enst2ensg
Mm.tnumber
Doelken, L., Ruzsics, Z., Raedle, B., Friedel, C. C., Zimmer, R., Mages, J., Hoffmann, R., Dickinson, P., Forster, T., Ghazal, P., & Koszinowski, U. H. (2008). High-resolution gene expression profiling for simultaneous kinetic parameter analysis of RNA synthesis and decay. RNA 14(9), 1959-1972. E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
DTA.dynamic.estimate uses an experiment, given by a phenotype matrix, data matrix and the number of uridines for each gene to estimate synthesis and decay rate of the genes.
DTA.dynamic.estimate(phenomat = NULL,datamat = NULL,tnumber = NULL,ccl = NULL,mRNAs = NULL,reliable = NULL,mediancenter = TRUE,usefractions = "LandT",LtoTratio = NULL,ratiomethod = "tls",largest = 5,weighted = TRUE,relevant = NULL,check = TRUE,error = TRUE,samplesize = 1000,confidence.range = c(0.025,0.975),bicor = TRUE,condition = "",upper = 700,lower = 500,save.plots = FALSE,resolution = 1,folder = NULL,fileformat = "jpeg",totaloverwt = 1,sr.vs.dr.folds.lims = c(-5,5),te.vs.to.folds.lims = c(-6,6),robust = FALSE,clusters = "sr",ranktime = NULL,upperquant = 0.8,lowerquant = 0.6,notinR = FALSE,RStudio = FALSE,simulation = FALSE,sim.object = NULL)
DTA.dynamic.estimate(phenomat = NULL,datamat = NULL,tnumber = NULL,ccl = NULL,mRNAs = NULL,reliable = NULL,mediancenter = TRUE,usefractions = "LandT",LtoTratio = NULL,ratiomethod = "tls",largest = 5,weighted = TRUE,relevant = NULL,check = TRUE,error = TRUE,samplesize = 1000,confidence.range = c(0.025,0.975),bicor = TRUE,condition = "",upper = 700,lower = 500,save.plots = FALSE,resolution = 1,folder = NULL,fileformat = "jpeg",totaloverwt = 1,sr.vs.dr.folds.lims = c(-5,5),te.vs.to.folds.lims = c(-6,6),robust = FALSE,clusters = "sr",ranktime = NULL,upperquant = 0.8,lowerquant = 0.6,notinR = FALSE,RStudio = FALSE,simulation = FALSE,sim.object = NULL)
phenomat |
A phenotype matrix, containing the design of the experiment as produced by |
datamat |
A matrix, containing the measurements from U, L and T, according to the design given in phenomat. Matrix should only contain the rows of phenomat as columns. |
tnumber |
Integer vector, containing the numbers of uridines. Elements should have the rownames of datamat. |
ccl |
The cell cycle length of the cells. |
mRNAs |
Estimated number of mRNAs in a cell (optional). |
reliable |
Vector of 'reliable' genes, which are used for parameter estimation. |
mediancenter |
Should the quotient Labeled/Total resp. Unlabeled/Total be rescaled to a common median over it's replicates before building the genewise median. |
usefractions |
From which fractions should the decay rate be calculated: "LandT", "UandT" or "both". |
LtoTratio |
Coefficient to rescale Labeled/Total. Is estimated from the data, if not specified. See ratiomethod. |
ratiomethod |
Choose the regression method to be used, possible methods are: "tls", "bias" and "lm". For details, see supplemental material of Sun et al. (see references). |
largest |
Percentage of largest residues from the first regression not to be used in the second regression step. For details, see supplemental material of Sun et al. (see references). |
weighted |
Should the regression be weighted with 1/(Total^2 + median(Total))? |
relevant |
Choose the arrays to be used for halflives calculation, vector due to nr (=replicate number) in phenomat. |
check |
If check = TRUE, control messages and plots will be generated. |
error |
If TRUE, the measurement error is assessed by means of an error model and resampling to gain confidence regions. |
samplesize |
Error model samplesize for resampling. |
confidence.range |
Confidence region for error model as quantiles. Interval should be between 0 and 1. |
bicor |
Should the labeling bias be corrected? |
condition |
String, to be added to the plotnames. |
upper |
Upper bound for labeling bias estimation. For details, see supplemental material of Sun et al. (see references). |
lower |
Lower bound for labeling bias estimation. For details, see supplemental material of Sun et al. (see references). |
save.plots |
If save.plots = TRUE, control plots will be saved. |
resolution |
Resolution scaling factor for plotting. |
folder |
Path to the folder, where to save the plots. |
fileformat |
Fileformat for plots to be saved. See |
totaloverwt |
Will be available in the very near future for comparative DTA data. |
sr.vs.dr.folds.lims |
Limits of the folds plot (dr vs sr). |
te.vs.to.folds.lims |
Limits of the folds plot (LT vs LE). |
robust |
If robust = TRUE, LE resp. LT is chosen instead of sr resp. dr. |
clusters |
should the dr vs sr folds be plotted with clusters, choose 'sr', 'dr' for cluster selection or 'none' to omit it |
ranktime |
at which time should the rankgain be calculated, default is the last column |
upperquant |
upper quantile for cluster selection |
lowerquant |
lower quantile for cluster selection |
notinR |
Should plots be not plotted in R. |
RStudio |
For RStudio users. Suppresses the opening of a new device, as RStudio allows only one. |
simulation |
True, if data was generated by |
sim.object |
Simulation object created by |
DTA.dynamic.estimate
returns a list, where each entry contains the estimation results for all replicates of one timecourse timepoint. Each result contains the following entries
triples |
Mapping of each fraction and experiment to its corresponding column in the data matrix. |
plabel |
The labeling efficiency. For details, see the vignette. |
LtoTratio |
Estimated ratio of labeled to total fraction. |
UtoTratio |
Estimated ratio of unlabeled to total fraction. |
LtoUratio |
Estimated ratio of labeled to unlabeled fraction. |
correcteddatamat |
Labeling bias corrected data matrix. |
drmat |
Decay rates for each replicate. The last column gives the median decay rates. |
dr |
Median decay rates. The last column of drmat. |
dr.confidence |
Confidence regions of decay rates. |
hlmat |
Half-lives for each replicate. The last column gives the median half-lifes. |
hl |
Median half-lives. The last column of hlmat. |
hl.confidence |
Confidence regions of half-lives. |
TEmat |
Total expression for each replicate. The last column gives the median total expression values. |
TE |
Median total expression values. The last column of TEmat. |
TE.confidence |
Confidence regions of total expression values. |
LEmat |
Labeled expression for each replicate. The last column gives the median labeled expression values. |
LE |
Median labeled expression values. The last column of LEmat. |
LE.confidence |
Confidence regions of labeled expression values. |
UEmat |
Unlabeled expression for each replicate. The last column gives the median unlabeled expression values. (Only if unlabeled values exist in the experiment) |
UE |
Median unlabeled expression values. The last column of UEmat. (Only if unlabeled values exist in the experiment) |
UE.confidence |
Confidence regions of unlabeled expression values. |
srmat |
Synthesis rates for each replicate. The last column gives the median synthesis rates. |
sr |
Median synthesis rates. The last column of srmat. |
sr.confidence |
Confidence regions of synthesis rates. |
LtoTmat |
Labeled to total ratio for each replicate. The last column gives the median labeled to total ratios. |
LtoT |
Median labeled to total ratios. The last column of LtoTmat. |
LtoT.confidence |
Confidence regions of labeled to total ratios. |
UtoTmat |
Unlabeled to total ratio for each replicate. The last column gives the median unlabeled to total ratios. |
UtoT |
Median unlabeled to total ratios. The last column of UtoTmat. |
UtoT.confidence |
Confidence regions of unlabeled to total ratios. |
Rsrmat |
Rescaled synthesis rates for each replicate, if parameter |
Rsr |
Rescaled median synthesis rates. The last column of Rsrmat. |
globaldrmat |
Decay rate for each replicate. Reciprocally weighted by the total expression. Last element contains (weighted) median decay rate. |
globaldr |
(Weighted) median decay rate. |
Bjoern Schwalb [email protected]
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011. M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review. B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics.
dataPath = system.file("data", package="DTA") load(file.path(dataPath, "Miller2011dynamic.RData")) ### for control plots set 'check = TRUE' ### res = DTA.dynamic.estimate(Sc.phenomat.dynamic,Sc.datamat.dynamic,Sc.tnumber,ccl = 150,mRNAs = 60000,reliable = Sc.reliable.dynamic,LtoTratio = rep(0.1,7),check = FALSE)
dataPath = system.file("data", package="DTA") load(file.path(dataPath, "Miller2011dynamic.RData")) ### for control plots set 'check = TRUE' ### res = DTA.dynamic.estimate(Sc.phenomat.dynamic,Sc.datamat.dynamic,Sc.tnumber,ccl = 150,mRNAs = 60000,reliable = Sc.reliable.dynamic,LtoTratio = rep(0.1,7),check = FALSE)
DTA.dynamic.generate produces the phenotype matrix and the matrix containing the simulated data according to the given parameters.
DTA.dynamic.generate(duration = 60,lab.duration = 6,tnumber = NULL,plabel = NULL,nrgenes = 5000,mediantime.halflives = 12,mediantime.synthesisrates = 18,n = 1,ccl = NULL,check = TRUE,plots = FALSE,save.plots = FALSE,folder = NULL,condition = "",addformat = NULL,sdnoise = 0.075,nobias = FALSE,unspecific.LtoU = 0,unspec.LtoU.weighted = FALSE,unspecific.UtoL = 0,unspec.UtoL.weighted = FALSE,mu.values.mat = NULL,mu.breaks.mat = NULL,lambda.values.mat = NULL,lambda.breaks.mat = NULL,truehalflives = NULL,truesynthesisrates = NULL,genenames = NULL)
DTA.dynamic.generate(duration = 60,lab.duration = 6,tnumber = NULL,plabel = NULL,nrgenes = 5000,mediantime.halflives = 12,mediantime.synthesisrates = 18,n = 1,ccl = NULL,check = TRUE,plots = FALSE,save.plots = FALSE,folder = NULL,condition = "",addformat = NULL,sdnoise = 0.075,nobias = FALSE,unspecific.LtoU = 0,unspec.LtoU.weighted = FALSE,unspecific.UtoL = 0,unspec.UtoL.weighted = FALSE,mu.values.mat = NULL,mu.breaks.mat = NULL,lambda.values.mat = NULL,lambda.breaks.mat = NULL,truehalflives = NULL,truesynthesisrates = NULL,genenames = NULL)
duration |
duration of the whole time course (min) |
lab.duration |
labeling duration for single experiments (min) |
tnumber |
Integer vector containing the number of uridine residues for each gene. If NULL, tnumber is sampled from an F-distribution within the function. |
plabel |
The labeling efficiency. If NULL, plabel is set to 0.005 within the function. For details, see supplemental material of Sun et al. (see references). |
nrgenes |
The number of genes the simulated experiment will have (will be cropped if it exceeds the length of tnumber). |
mediantime.halflives |
the median of the half life distribution |
mediantime.synthesisrates |
the median of the synthesis rates distribution (counts/cell/cellcycle) |
n |
the number of cells N(0) |
ccl |
The cell cycle length (in minutes). |
check |
if check=TRUE, control messages will be generated |
plots |
if plots = TRUE, control plots will be plotted |
save.plots |
if save.plots = TRUE, control plots will be saved |
folder |
folder, where to save the plots |
condition |
to be added to the plotnames |
addformat |
additional fileformat for plots to be saved |
sdnoise |
The amount of measurement noise (proportional to expression strength). |
nobias |
Should a labeling bias be added? |
unspecific.LtoU |
Proportion of labeled RNAs that unspecifically end up in the unlabeled fraction. |
unspec.LtoU.weighted |
Should unspecific proportion of labeled to unlabeled depend linearly on the length of the RNA? |
unspecific.UtoL |
Proportion of unlabeled RNAs that unspecifically end up in the labeled fraction. |
unspec.UtoL.weighted |
Should unspecific proportion of unlabeled to labeled depend linearly on the length of the RNA? |
mu.values.mat |
if the data should be generated using given synthesis rates, this matrix must contain the respective values for each gene |
mu.breaks.mat |
timepoints of synthesis rate changes, this matrix must contain the respective values for each gene, only needed when mu.values.mat is given (one column less than mu.values.mat) |
lambda.values.mat |
if the data should be generated using given decay rates, this matrix must contain the respective values for each gene |
lambda.breaks.mat |
timepoints of decay rate changes, this matrix must contain the respective values for each gene, only needed when lambda.values.mat is given (one column less than lambda.values.mat) |
truehalflives |
If the data should be generated using a given half-life distribution, this vector must contain the respective values for each gene. |
truesynthesisrates |
If the data should be generated using a given synthesis rates distribution, this vector must contain the respective values for each gene |
genenames |
An optional list of gene names. |
DTA.dynamic.generate returns a list, containing the following entries
phenomat |
A matrix, containing the design of the experiment as produced by |
datamat |
A matrix, containing the simulated measurements from U, L and T, according to the design given in phenomat. |
tnumber |
Integer vector containing the number of uridine residues for each gene. |
ccl |
The cell cycle length (in minutes). |
truemus |
A vector, containing the true synthesis rates. |
truemusaveraged |
A vector, containing the true synthesis rates, averaged over the labeling period. |
truelambdas |
A vector, containing the true decay rates. |
truelambdasaveraged |
A vector, containing the true decay rates, averaged over the labeling period. |
truehalflives |
A vector, containing the true half-lives. |
truehalflivesaveraged |
A vector, containing the true half-lives, averaged over the labeling period. |
trueplabel |
The true labeling efficiency. For details, see supplemental material of Sun et al. (see references). |
truecomplete |
A vector, containing the true amount of total RNA. |
truelambdas |
A vector, containing the true decay rates. |
truemus |
A vector, containing the true synthesis rates. |
truehalflives |
A vector, containing the true half-lives. |
trueplabel |
The true labeling efficiency. For details, see supplemental material of Miller et al. (see references). |
truear |
The true parameter ar. For details, see supplemental material of Miller et al. (see references). |
truebr |
The true parameter br. For details, see supplemental material of Miller et al. (see references). |
truecr |
The true parameter cr. For details, see supplemental material of Miller et al. (see references). |
truecrbyar |
The true parameter cr/ar. For details, see supplemental material of Miller et al. (see references). |
truecrbybr |
The true parameter cr/br. For details, see supplemental material of Miller et al. (see references). |
truebrbyar |
The true parameter br/ar. For details, see supplemental material of Miller et al. (see references). |
trueLasymptote |
The true parameter asymptote (labeled bias). For details, see supplemental material of Miller et al. (see references). |
trueUasymptote |
The true parameter asymptote (unlabeled bias). For details, see supplemental material of Miller et al. (see references). |
Bjoern Schwalb [email protected]
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Dolken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011. M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review. B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics.
nrgenes = 5000 truesynthesisrates = rf(nrgenes,5,5)*18 steady = rep(1,nrgenes) shock = 1/pmax(rnorm(nrgenes,mean = 8,sd = 4),1) induction = pmax(rnorm(nrgenes,mean = 8,sd = 4),1) changes.mat = cbind(steady,shock,shock*induction) mu.values.mat = changes.mat*truesynthesisrates mu.breaks.mat = cbind(rep(12,nrgenes),rep(18,nrgenes)) truehalflives = rf(nrgenes,15,15)*12 truelambdas = log(2)/truehalflives changes.mat = cbind(steady,shock,shock*induction,steady) lambda.values.mat = changes.mat*truelambdas lambda.breaks.mat = cbind(rep(12,nrgenes),rep(18,nrgenes),rep(27,nrgenes)) ### it takes several min to build sim.object (depends on the number of genes 'nrgenes') ### sim.object = DTA.dynamic.generate(duration = 36,lab.duration = 6,nrgenes = nrgenes,mu.values.mat = mu.values.mat,mu.breaks.mat = mu.breaks.mat,lambda.values.mat = lambda.values.mat,lambda.breaks.mat = lambda.breaks.mat) ### for control plots set 'check = TRUE' ### res = DTA.dynamic.estimate(simulation = TRUE,sim.object = sim.object,check = FALSE)
nrgenes = 5000 truesynthesisrates = rf(nrgenes,5,5)*18 steady = rep(1,nrgenes) shock = 1/pmax(rnorm(nrgenes,mean = 8,sd = 4),1) induction = pmax(rnorm(nrgenes,mean = 8,sd = 4),1) changes.mat = cbind(steady,shock,shock*induction) mu.values.mat = changes.mat*truesynthesisrates mu.breaks.mat = cbind(rep(12,nrgenes),rep(18,nrgenes)) truehalflives = rf(nrgenes,15,15)*12 truelambdas = log(2)/truehalflives changes.mat = cbind(steady,shock,shock*induction,steady) lambda.values.mat = changes.mat*truelambdas lambda.breaks.mat = cbind(rep(12,nrgenes),rep(18,nrgenes),rep(27,nrgenes)) ### it takes several min to build sim.object (depends on the number of genes 'nrgenes') ### sim.object = DTA.dynamic.generate(duration = 36,lab.duration = 6,nrgenes = nrgenes,mu.values.mat = mu.values.mat,mu.breaks.mat = mu.breaks.mat,lambda.values.mat = lambda.values.mat,lambda.breaks.mat = lambda.breaks.mat) ### for control plots set 'check = TRUE' ### res = DTA.dynamic.estimate(simulation = TRUE,sim.object = sim.object,check = FALSE)
DTA.estimate uses an experiment, given by a phenotype matrix, data matrix and the number of uridines for each gene to estimate synthesis and decay rate of the genes.
DTA.estimate(phenomat = NULL,datamat = NULL,tnumber = NULL,reliable = NULL,ccl = NULL,mRNAs = NULL,mediancenter = TRUE,usefractions = "LandT",LtoTratio = NULL,ratiomethod = "tls",largest = 5,weighted = TRUE,relevant = NULL,upper = 700,lower = 500,error = TRUE,samplesize = 1000,confidence.range = c(0.025,0.975),bicor = TRUE,check = TRUE,condition = "",save.plots = FALSE,resolution = 1,notinR = FALSE,RStudio = FALSE,folder = NULL,fileformat = "jpeg",totaloverwt = 1,simulation = FALSE,sim.object = NULL)
DTA.estimate(phenomat = NULL,datamat = NULL,tnumber = NULL,reliable = NULL,ccl = NULL,mRNAs = NULL,mediancenter = TRUE,usefractions = "LandT",LtoTratio = NULL,ratiomethod = "tls",largest = 5,weighted = TRUE,relevant = NULL,upper = 700,lower = 500,error = TRUE,samplesize = 1000,confidence.range = c(0.025,0.975),bicor = TRUE,check = TRUE,condition = "",save.plots = FALSE,resolution = 1,notinR = FALSE,RStudio = FALSE,folder = NULL,fileformat = "jpeg",totaloverwt = 1,simulation = FALSE,sim.object = NULL)
phenomat |
A phenotype matrix, containing the design of the experiment as produced by |
datamat |
A matrix, containing the measurements from U, L and T, according to the design given in phenomat. Matrix should only contain the rows of phenomat as columns. |
tnumber |
Integer vector, containing the numbers of uridine residues for each transcript. Elements should have the rownames of datamat. |
ccl |
The cell cycle length of the cells (optional). Is not modeled, if not set. |
mRNAs |
Estimated number of mRNAs in a cell (optional). |
reliable |
Vector of 'reliable' genes, which are used for parameter estimation. |
mediancenter |
Should the quotient Labeled/Total resp. Unlabeled/Total be rescaled to a common median over it's replicates before building the genewise median. |
usefractions |
From which fractions should the decay rate be calculated: "LandT", "UandT" or "both". |
LtoTratio |
Coefficient to rescale Labeled/Total. Is estimated from the data, if not specified. See ratiomethod. Altering this parameter leads to a altered median half-life. For details, see supplemental material of Sun et al. (see references). |
ratiomethod |
Choose the regression method to be used, possible methods are: "tls", "bias" and "lm". For details, see supplemental material of Sun et al. (see references). Method to estimate the parameter |
largest |
Percentage of largest residues from the first regression not to be used in the second regression step. For details, see supplemental material of Sun et al. (see references). |
weighted |
Should the regression be weighted with 1/(Total^2 + median(Total))? |
relevant |
Choose the arrays to be used for halflives calculation, vector due to nr (=replicate number) in phenomat. If not set, all arrays are used. |
check |
If check = TRUE, control messages and plots will be generated. |
error |
If TRUE, the measurement error is assessed by means of an error model and resampling to gain confidence regions. |
samplesize |
Error model samplesize for resampling. |
confidence.range |
Confidence region for error model as quantiles. Interval should be between 0 and 1. |
bicor |
Should the labeling bias be corrected? |
condition |
String, to be added to the plotnames if saved. |
upper |
Upper bound for labeling bias estimation. For details, see supplemental material of Sun et al. (see references). |
lower |
Lower bound for labeling bias estimation. For details, see supplemental material of Sun et al. (see references). |
save.plots |
If save.plots = TRUE, control plots will be saved. Please check folder writability. |
resolution |
Resolution scaling factor for plotting. (Scaled with 72dpi.) |
notinR |
If TRUE, plots are not plotted in R. |
RStudio |
For RStudio users. Suppresses the opening of a new device, as RStudio allows only one. |
folder |
Path to the folder, where to save the plots. Needs to be writable. |
fileformat |
Fileformat for plots to be saved. See |
totaloverwt |
Only needed when |
simulation |
True, if data was generated by |
sim.object |
Simulation object created by |
DTA.estimate
returns a list, where each entry contains the estimation results for all replicates of one labeling time. Each result contains the following entries
triples |
Mapping of each fraction and experiment to its corresponding column in the data matrix. |
plabel |
The labeling efficiency. For details, see supplemental material of Sun et al. (see references). |
LtoTratio |
Estimated ratio of labeled to total fraction. |
UtoTratio |
Estimated ratio of unlabeled to total fraction. |
LtoUratio |
Estimated ratio of labeled to unlabeled fraction. |
correcteddatamat |
Labeling bias corrected data matrix. |
drmat |
Decay rates for each replicate. The last column gives the median decay rates. |
dr |
Median decay rates. The last column of drmat. |
dr.confidence |
Confidence regions of decay rates. |
hlmat |
Half-lives for each replicate. The last column gives the median half-lifes. |
hl |
Median half-lives. The last column of hlmat. |
hl.confidence |
Confidence regions of half-lives. |
TEmat |
Total expression for each replicate. The last column gives the median total expression values. |
TE |
Median total expression values. The last column of TEmat. |
TE.confidence |
Confidence regions of total expression values. |
LEmat |
Labeled expression for each replicate. The last column gives the median labeled expression values. |
LE |
Median labeled expression values. The last column of LEmat. |
LE.confidence |
Confidence regions of labeled expression values. |
UEmat |
Unlabeled expression for each replicate. The last column gives the median unlabeled expression values. (Only if unlabeled values exist in the experiment) |
UE |
Median unlabeled expression values. The last column of UEmat. (Only if unlabeled values exist in the experiment) |
UE.confidence |
Confidence regions of unlabeled expression values. |
srmat |
Synthesis rates for each replicate. The last column gives the median synthesis rates. |
sr |
Median synthesis rates. The last column of srmat. |
sr.confidence |
Confidence regions of synthesis rates. |
LtoTmat |
Labeled to total ratio for each replicate. The last column gives the median labeled to total ratios. |
LtoT |
Median labeled to total ratios. The last column of LtoTmat. |
LtoT.confidence |
Confidence regions of labeled to total ratios. |
UtoTmat |
Unlabeled to total ratio for each replicate. The last column gives the median unlabeled to total ratios. |
UtoT |
Median unlabeled to total ratios. The last column of UtoTmat. |
UtoT.confidence |
Confidence regions of unlabeled to total ratios. |
Rsrmat |
Rescaled synthesis rates for each replicate, if parameter |
Rsr |
Rescaled median synthesis rates. The last column of Rsrmat. |
globaldrmat |
Decay rate for each replicate. Reciprocally weighted by the total expression. Last element contains (weighted) median decay rate. |
globaldr |
(Weighted) median decay rate. |
Bjoern Schwalb [email protected]
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011. M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review. B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics.
dataPath = system.file("data", package="DTA") load(file.path(dataPath, "Miller2011.RData")) ### for control plots set 'check = TRUE' ### res = DTA.estimate(Sc.phenomat,Sc.datamat,Sc.tnumber,ccl = 150,mRNAs = 60000,reliable = Sc.reliable,check = FALSE)
dataPath = system.file("data", package="DTA") load(file.path(dataPath, "Miller2011.RData")) ### for control plots set 'check = TRUE' ### res = DTA.estimate(Sc.phenomat,Sc.datamat,Sc.tnumber,ccl = 150,mRNAs = 60000,reliable = Sc.reliable,check = FALSE)
DTA.generate produces the phenotype matrix and the matrix containing the simulated data according to the given parameters.
DTA.generate(timepoints, tnumber = NULL, plabel = NULL, nrgenes = 5000, mediantime = 12, ccl = 150, delaytime = 0, sdnoise = 0.075, nobias = FALSE, unspecific.LtoU = 0, unspec.LtoU.weighted = FALSE, unspecific.UtoL = 0, unspec.UtoL.weighted = FALSE, truehalflives = NULL, truecomplete = NULL, genenames = NULL, cDTA = FALSE)
DTA.generate(timepoints, tnumber = NULL, plabel = NULL, nrgenes = 5000, mediantime = 12, ccl = 150, delaytime = 0, sdnoise = 0.075, nobias = FALSE, unspecific.LtoU = 0, unspec.LtoU.weighted = FALSE, unspecific.UtoL = 0, unspec.UtoL.weighted = FALSE, truehalflives = NULL, truecomplete = NULL, genenames = NULL, cDTA = FALSE)
timepoints |
Integer vector containing the labeling times for which the samples should be generated. |
tnumber |
Integer vector containing the number of uridine residues for each gene. If NULL, tnumber is sampled from an F-distribution within the function. |
plabel |
The labeling efficiency. If NULL, plabel is set to 0.005 within the function. For details, see supplemental material of Sun et al. (see references). |
nrgenes |
The number of genes the simulated experiment will have (will be cropped if it exceeds the length of tnumber). |
mediantime |
The median of the randomly drawn half-life distribution. |
ccl |
The cell cycle length (in minutes). |
delaytime |
Estimates the delay between the moment of 4sU/4tU labeling and actual incorporation of it into mRNA. |
sdnoise |
The amount of measurement noise (proportional to expression strength). |
nobias |
Should a labeling bias be added? |
unspecific.LtoU |
Proportion of labeled RNAs that unspecifically end up in the unlabeled fraction. |
unspec.LtoU.weighted |
Should unspecific proportion of labeled to unlabeled depend linearly on the length of the RNA? |
unspecific.UtoL |
Proportion of unlabeled RNAs that unspecifically end up in the labeled fraction. |
unspec.UtoL.weighted |
Should unspecific proportion of unlabeled to labeled depend linearly on the length of the RNA? |
truehalflives |
If the data should be generated using a given half-life distribution, this vector must contain the respective values for each gene. |
truecomplete |
If the data should be generated using a given expression distribution, this vector must contain the respective values for each gene. |
genenames |
An optional list of gene names. |
cDTA |
cDTA = FALSE does not rescale L and U. |
DTA.generate returns a list, containing the following entries
phenomat |
A matrix, containing the design of the experiment as produced by |
datamat |
A matrix, containing the simulated measurements from U, L and T, according to the design given in phenomat. |
tnumber |
Integer vector containing the number of uridine residues for each gene. |
ccl |
The cell cycle length (in minutes). |
truecomplete |
A vector, containing the true amount of total RNA. |
truelambdas |
A vector, containing the true decay rates. |
truemus |
A vector, containing the true synthesis rates. |
truehalflives |
A vector, containing the true half-lives. |
trueplabel |
The true labeling efficiency. For details, see supplemental material of Miller et al. (see references). |
truear |
The true parameter ar. For details, see supplemental material of Miller et al. (see references). |
truebr |
The true parameter br. For details, see supplemental material of Miller et al. (see references). |
truecr |
The true parameter cr. For details, see supplemental material of Miller et al. (see references). |
truecrbyar |
The true parameter cr/ar. For details, see supplemental material of Miller et al. (see references). |
truecrbybr |
The true parameter cr/br. For details, see supplemental material of Miller et al. (see references). |
truebrbyar |
The true parameter br/ar. For details, see supplemental material of Miller et al. (see references). |
trueLasymptote |
The true parameter asymptote (labeled bias). For details, see supplemental material of Miller et al. (see references). |
trueUasymptote |
The true parameter asymptote (unlabeled bias). For details, see supplemental material of Miller et al. (see references). |
Bjoern Schwalb [email protected]
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011. M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review. B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics.
sim.object = DTA.generate(timepoints=rep(c(6,12),2)) ### for control plots set 'check = TRUE' ### res.sim = DTA.estimate(ratiomethod = "bias",simulation = TRUE,sim.object = sim.object,check = FALSE)
sim.object = DTA.generate(timepoints=rep(c(6,12),2)) ### for control plots set 'check = TRUE' ### res.sim = DTA.estimate(ratiomethod = "bias",simulation = TRUE,sim.object = sim.object,check = FALSE)
DTA.map.it
can map different kinds of identifiers in a matrix or a vector given by mapping vector.
DTA.map.it(mat,map = NULL,check = TRUE)
DTA.map.it(mat,map = NULL,check = TRUE)
mat |
Matrix or vector with numerical entries. |
map |
Vector of identifiers to map to, named by identifiers to map from. |
check |
Should check protocol be printed. |
Bjoern Schwalb [email protected]
M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review. B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics.
### see vignette examples or reference: ### B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. ### Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics, in revision.
### see vignette examples or reference: ### B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. ### Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics, in revision.
DTA.normalize
can normalize expression values from a certain species to the median of values from a reference species.
DTA.normalize(mat,reliable = NULL,logscale = FALSE,protocol = FALSE,center = FALSE)
DTA.normalize(mat,reliable = NULL,logscale = FALSE,protocol = FALSE,center = FALSE)
mat |
Expression matrix. |
reliable |
The rows to be used, i.e. identifiers of the reference species to normalize on. |
logscale |
Is the matrix in log-scale ? |
protocol |
Should a protocol be printed ? |
center |
Should the center be 0 (log-scale) or 1 (absolute scale). Otherwise the median of the medians is taken. |
Bjoern Schwalb [email protected]
M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review. B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics.
### see vignette examples or reference: ### B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. ### Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics, in revision.
### see vignette examples or reference: ### B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. ### Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics, in revision.
DTA.phenomat
creates a phenomat for a given experimental design, i.e. used labeling times.
DTA.phenomat(timepoints, timecourse = NULL)
DTA.phenomat(timepoints, timecourse = NULL)
timepoints |
The respective labeling times of the measured samples. |
timecourse |
Vector giving the order for timecourse DTA data. |
A matrix, containing the design of the experiment. Columns are name, fraction (U=unlabebeld, L=labeled, T=total), time and nr (=replicate number). Rows represent individual experiments. For timecourse data, an additional column of the order of the underlying timecourse data can be added via timecourse
.
Bjoern Schwalb [email protected]
### phenomat for 2 replicates of 6 and 12 min labeling duration resp. DTA.phenomat(c(6,12)) ### phenomat for three adjacent timepoints measured in 2 replicates DTA.phenomat(rep(6,6),timecourse = 1:3)
### phenomat for 2 replicates of 6 and 12 min labeling duration resp. DTA.phenomat(c(6,12)) ### phenomat for three adjacent timepoints measured in 2 replicates DTA.phenomat(rep(6,6),timecourse = 1:3)
DTA.plot.it
can save plots in any format and any quality in addition to show them in R devices
DTA.plot.it(filename,sw = 1,sh = 1,sres = 1,plotsfkt,ww = 7,wh = 7,pointsize = 12,dev.pointsize = 8,paper = "special",quality = 100,units = "px",bg = "white",fileformat = "jpeg",saveit = FALSE,notinR = FALSE,RStudio = FALSE,addformat = NULL)
DTA.plot.it(filename,sw = 1,sh = 1,sres = 1,plotsfkt,ww = 7,wh = 7,pointsize = 12,dev.pointsize = 8,paper = "special",quality = 100,units = "px",bg = "white",fileformat = "jpeg",saveit = FALSE,notinR = FALSE,RStudio = FALSE,addformat = NULL)
filename |
Name of the plot to be saved without the format type suffix. |
sw |
Scaling factor of width. Scaled with 480px. |
sh |
Scaling factor of height. Scaled with 480px. |
sres |
Scaling factor of the resolution. Scaled with 72dpi. |
plotsfkt |
Function of plots to be plotted. |
ww |
Width of window. Needed only for plotting in R or if filformat = "pdf" or "ps". See pdf or ps. |
wh |
Height of window. Needed only for plotting in R or if filformat = "pdf" or "ps". See pdf or ps. |
pointsize |
The default pointsize of plotted text, interpreted as big points (1/72 inch) for plots to be saved. |
dev.pointsize |
Pointsize of plotted text, interpreted as big points (1/72 inch) for display in R. |
paper |
Needed only if filformat = "pdf" or "ps". See pdf or ps. |
quality |
Needed only if filformat = "jpeg". See jpeg. |
units |
Needed only if filformat = "jpeg", "png", "bmp" or "tiff". See corresponding function. |
bg |
Backgroundcolor. |
fileformat |
Save the plot as "jpeg", "png", "bmp", "tiff", "ps" or "pdf". |
saveit |
Should plot be saved. |
notinR |
Should plot be not plotted in R. |
RStudio |
For RStudio users. Suppresses the opening of a new device, as RStudio allows only one. |
addformat |
Should plot be saved additionally in another format, "jpeg", "png", "bmp", "tiff", "ps" or "pdf". |
Bjoern Schwalb [email protected]
plotsfkt = function(){ par(mfrow = c(1,2)) plot(1:10) plot(10:1) } DTA.plot.it(filename = "test",plotsfkt = plotsfkt,saveit = TRUE) dev.off()
plotsfkt = function(){ par(mfrow = c(1,2)) plot(1:10) plot(10:1) } DTA.plot.it(filename = "test",plotsfkt = plotsfkt,saveit = TRUE) dev.off()
This matrix contains the RNA intensity values for each gene across each RNA fraction and their replicate measurements of the Homo Sapiens DTA experiment from Doelken et al.
Hs.datamat
Hs.datamat
The column names of the matrix give the cel-file name and the row names the Ensembl gene IDs.
Doelken, L., Ruzsics, Z., Raedle, B., Friedel, C. C., Zimmer, R., Mages, J., Hoffmann, R., Dickinson, P., Forster, T., Ghazal, P., & Koszinowski, U. H. (2008). High-resolution gene expression profiling for simultaneous kinetic parameter analysis of RNA synthesis and decay. RNA 14(9), 1959-1972.
Mapping from Ensembl transcript IDs to Ensembl gene IDs of Homo Sapiens.
Hs.enst2ensg
Hs.enst2ensg
Vector gives the Ensembl gene IDs, names the Ensembl transcript IDs.
E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
The phenotype matrix Hs.phenomat
contains information about the experimental design. It is comprised of the filename, the type of RNA fraction measured (T, U or L), the labeling time and the replicate number.
Hs.phenomat
Hs.phenomat
The phenomat is a matrix comprised of the file name, the type of RNA fraction mesasured (T, U or L, fraction
column), the labeling time (time
,timeframe
column) and the replicate number (nr
column). Rows in this matrix represent the individual experiments.
Doelken, L., Ruzsics, Z., Raedle, B., Friedel, C. C., Zimmer, R., Mages, J., Hoffmann, R., Dickinson, P., Forster, T., Ghazal, P., & Koszinowski, U. H. (2008). High-resolution gene expression profiling for simultaneous kinetic parameter analysis of RNA synthesis and decay. RNA 14(9), 1959-1972.
Ensembl gene IDs, that passed certain criteria among the Homo Sapiens Doelken et al. DTA experiment to be considered valid for parameter estimation. For details, see vignette.
Hs.reliable
Hs.reliable
Vector of Ensembl gene IDs that can be passed to DTA.estimate
for parameter estimation.
Doelken, L., Ruzsics, Z., Raedle, B., Friedel, C. C., Zimmer, R., Mages, J., Hoffmann, R., Dickinson, P., Forster, T., Ghazal, P., & Koszinowski, U. H. (2008). High-resolution gene expression profiling for simultaneous kinetic parameter analysis of RNA synthesis and decay. RNA 14(9), 1959-1972.
The amount of thymines in the cDNA of each transcript of all Homo Sapiens Ensembl transcript IDs, to assess the uridine-dependent labeling bias and eventually correct for it.
Hs.tnumber
Hs.tnumber
Vector gives the number of thymines in the cDNA (uridine residues in RNA) of each Ensembl transcript ID.
E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
R object contains all relevant *.RData files needed for the DTA.estimate
function. For example, see vignette.
Miller2011
Miller2011
R object contains the following *.RData files:
Sc.phenomat
Sc.datamat
Sc.reliable
Sc.tnumber
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011. E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
R object contains all relevant *.RData files needed for the DTA.estimate
function. For example, see vignette.
Miller2011dynamic
Miller2011dynamic
R object contains the following *.RData files:
Sc.phenomat.dynamic
Sc.datamat.dynamic
Sc.reliable.dynamic
Sc.tnumber
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011. E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
This matrix contains the RNA intensity values for each gene across each RNA fraction and their replicate measurements of the Mus Musculus DTA experiment from Doelken et al.
Mm.datamat
Mm.datamat
The column names of the matrix give the cel-file name and the row names the Ensembl gene IDs.
Doelken, L., Ruzsics, Z., Raedle, B., Friedel, C. C., Zimmer, R., Mages, J., Hoffmann, R., Dickinson, P., Forster, T., Ghazal, P., & Koszinowski, U. H. (2008). High-resolution gene expression profiling for simultaneous kinetic parameter analysis of RNA synthesis and decay. RNA 14(9), 1959-1972.
Mapping from Ensembl transcript IDs to Ensembl gene IDs of Mus Musculus.
Mm.enst2ensg
Mm.enst2ensg
Vector gives the Ensembl gene IDs, names the Ensembl transcript IDs.
E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
The phenotype matrix Mm.phenomat
contains information about the experimental design. It is comprised of the filename, the type of RNA fraction measured (T, U or L), the labeling time and the replicate number.
Mm.phenomat
Mm.phenomat
The phenomat is a matrix comprised of the file name, the type of RNA fraction mesasured (T, U or L, fraction
column), the labeling time (time
,timeframe
column) and the replicate number (nr
column). Rows in this matrix represent the individual experiments.
Doelken, L., Ruzsics, Z., Raedle, B., Friedel, C. C., Zimmer, R., Mages, J., Hoffmann, R., Dickinson, P., Forster, T., Ghazal, P., & Koszinowski, U. H. (2008). High-resolution gene expression profiling for simultaneous kinetic parameter analysis of RNA synthesis and decay. RNA 14(9), 1959-1972.
Ensembl gene IDs, that passed certain criteria among the Mus Musculus Doelken et al. DTA experiment to be considered valid for parameter estimation. For details, see vignette.
Mm.reliable
Mm.reliable
Vector of Ensembl gene IDs that can be passed to DTA.estimate
for parameter estimation.
Doelken, L., Ruzsics, Z., Raedle, B., Friedel, C. C., Zimmer, R., Mages, J., Hoffmann, R., Dickinson, P., Forster, T., Ghazal, P., & Koszinowski, U. H. (2008). High-resolution gene expression profiling for simultaneous kinetic parameter analysis of RNA synthesis and decay. RNA 14(9), 1959-1972.
The amount of thymines in the cDNA of each transcript of all Mus Musculus Ensembl transcript IDs, to assess the uridine-dependent labeling bias and eventually correct for it.
Mm.tnumber
Mm.tnumber
Vector gives the number of thymines in the cDNA (uridine residues in RNA) of each Ensembl transcript ID.
E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
The phenotype matrix Pol.phenomat
contains information about the experimental design. It is comprised of the filename, the type of RNA fraction measured (T, U or L), the labeling time and the replicate number.
Pol.phenomat
Pol.phenomat
The phenomat is a matrix comprised of the file name, the type of RNA fraction mesasured (T, U or L, fraction
column), the labeling time (time
,timeframe
column) and the replicate number (nr
column). Rows in this matrix represent the individual experiments.
M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review.
This matrix contains the RNA intensity values for each gene across each RNA fraction and their replicate measurements of the Saccharomyces Cerevisiae rpb1-N488D (Slow Polymerase) and wild-type cDTA experiment from Sun et al.
Raw.datamat
Raw.datamat
The column names of the matrix give the cel-file name and the row names the affymetrix IDs.
M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review.
Mapping from Affymetrix Yeast 2.0 IDs to Ensembl gene IDs of SaccharomycesCerevisiae.
Sc.affy2ensg
Sc.affy2ensg
Vector gives the Ensembl gene IDs, names the Affymetrix Yeast 2.0 IDs.
M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review.
This matrix contains the RNA intensity values for each gene across each RNA fraction and their replicate measurements of the Saccharomyces Cerevisiae wild-type DTA experiment from Miller et al.
Sc.datamat
Sc.datamat
The column names of the matrix give the cel-file name and the row names the Ensembl gene IDs.
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011.
This matrix contains the RNA intensity values for each gene across each RNA fraction and their replicate measurements of the Saccharomyces Cerevisiae salt stress DTA experiment from Miller et al.
Sc.datamat.dynamic
Sc.datamat.dynamic
The column names of the matrix give the cel-file name and the row names the Ensembl gene IDs.
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011.
Ensembl gene IDs, that passed certain criteria among the Saccharomyces Cerevisiae Sun et al. cDTA experiment to be considered valid for parameter estimation. For details, see Sun et al (Materials and Methods).
Sc.ensg.reliable
Sc.ensg.reliable
Vector of Ensembl gene IDs that can be passed to DTA.estimate
for parameter estimation.
M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review.
The phenotype matrix Sc.phenomat
contains information about the experimental design. It is comprised of the filename, the type of RNA fraction measured (T, U or L), the labeling time and the replicate number.
Sc.phenomat
Sc.phenomat
The phenomat is a matrix comprised of the file name, the type of RNA fraction mesasured (T, U or L, fraction
column), the labeling time (time
,timeframe
column) and the replicate number (nr
column). Rows in this matrix represent the individual experiments.
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011.
The phenotype matrix Sc.phenomat.dynamic
contains information about the experimental design. It is comprised of the filename, the type of RNA fraction measured (T, U or L), the labeling time, the replicate number and an additional number indicating the timecourse order.
Sc.phenomat.dynamic
Sc.phenomat.dynamic
The phenomat is a matrix comprised of the file name, the type of RNA fraction mesasured (T, U or L, fraction
column), the labeling time (time
,timeframe
column), the replicate number (nr
column) and a number indicating the timecourse order (timecourse
column). Rows in this matrix represent the individual experiments.
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011.
Ensembl gene IDs, that passed certain criteria among the Saccharomyces Cerevisiae Miller et al. wild-type DTA experiment to be considered valid for parameter estimation. For details, see supplemental material Miller et al.
Sc.reliable
Sc.reliable
Vector of Ensembl gene IDs that can be passed to DTA.estimate
for parameter estimation.
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011.
Ensembl gene IDs, that passed certain criteria among the Saccharomyces Cerevisiae Miller et al. salt stress DTA experiment to be considered valid for parameter estimation. For details, see supplemental material Miller et al.
Sc.reliable.dynamic
Sc.reliable.dynamic
Vector of Ensembl gene IDs that can be passed to DTA.estimate
for parameter estimation.
C. Miller, B. Schwalb, K. Maier, D. Schulz, S. Duemcke, B. Zacher, A. Mayer, J. Sydow, L. Marcinowski, L. Doelken, D. E. Martin, A. Tresch, and P. Cramer. Dynamic transcriptome analysis measures rates of mRNA synthesis and decay in yeast. Mol Syst Biol, 7:458, 2011.
ORF identifiers (Ensembl Gene ID) found to be associated with ribosome biogenesis, rRNA processing etc.
Sc.ribig.ensg
Sc.ribig.ensg
Vector of ORF identifiers (Ensembl Gene ID).
P. Jorgensen, I. RupeAi, J. R. Sharom, L. Schneper, J. R. Broach, and M. Tyers. A dynamic transcriptional network communicates growth potential to ribosome synthesis and critical cell size. Genes & Development, 18(20):2491-2505, October 2004.
ORF identifiers (Ensembl Gene ID) encoding for ribosomal protein genes.
Sc.rpg.ensg
Sc.rpg.ensg
Vector of ORF identifiers (Ensembl Gene ID).
A. Nakao, M. Yoshihama, and N. Kenmochi. RPG: the Ribosomal Protein Gene database. Nucleic acids research, 32(Database issue), January 2004.
ORF identifiers (Ensembl Gene ID) found to be associated with stress response by the iterative signature algorithm.
Sc.stress.ensg
Sc.stress.ensg
Vector of ORF identifiers (Ensembl Gene ID).
J. Ihmels, G. Friedlander, S. Bergmann, O. Sarig, Y. Ziv, and N. Barkai. Revealing modular organization in the yeast transcriptional network. Nature genetics, 31(4):370-377, August 2002.
ORF identifiers (Ensembl Gene ID) encoding for transcription factors.
Sc.tf.ensg
Sc.tf.ensg
Vector of ORF identifiers (Ensembl Gene ID).
K. D. MacIsaac, T. Wang, D. B. Gordon, D. K. Gifford, G. D. Stormo, and E. Fraenkel. An improved map of conserved regulatory sites for saccharomyces cerevisiae. BMC Bioinformatics, 7:113, 2006.
The amount of thymines in the cDNA of each transcript of all Saccharomyces Cerevisiae Ensembl transcript IDs (ORF identifier), to assess the uridine-dependent labeling bias and eventually correct for it.
Sc.tnumber
Sc.tnumber
Vector gives the number of thymines in the cDNA (uridine residues in RNA) of each Ensembl transcript ID.
E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
Ensembl gene IDs, that passed certain criteria among the Saccharomyces Cerevisiae Sun et al. cDTA experiment to be considered valid for cDTA normalization. For details, see Sun et al (Materials and Methods).
Sp.affy.reliable
Sp.affy.reliable
Vector of Schizosaccharomyces Pombe affymetrix IDs that can be passed to DTA.normalize
for cDTA normalization of the Saccharomyces Cerevisiae identifiers.
M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review.
The amount of thymines in the cDNA of each transcript of all Schizosaccharomyces Pombe Ensembl transcript IDs (ORF identifier), to assess the uridine-dependent labeling bias and eventually correct for it.
Sp.tnumber
Sp.tnumber
Vector gives the number of thymines in the cDNA (uridine residues in RNA) of each Ensembl transcript ID.
E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
R object contains all relevant *.RData files needed for the DTA.estimate
function. For example, see Schwalb et al.
Sun2011
Sun2011
R object contains the following *.RData files:
Raw.datamat
Sp.affy.reliable
Sc.affy2ensg
Wt.phenomat
Pol.phenomat
Sc.ensg.reliable
Sc.tnumber
M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review. B. Schwalb, B. Zacher, S. Duemcke, D. Martin, P. Cramer, A. Tresch. Measurement of genome-wide RNA synthesis and decay rates with Dynamic Transcriptome Analysis (DTA/cDTA). Bioinformatics. E. Birney, D. Andrews, M. Caccamo, Y. Chen, L. Clarke, G. Coates, T. Cox, F. Cunningham, V. Curwen, T. Cutts, T. Down, R. Durbin, X. M. Fernandez-Suarez, P. Flicek, S. Graef, M. Hammond, J. Herrero, K. Howe, V. Iyer, K. Jekosch, A. Kaehaeri, A. Kasprzyk, D. Keefe, F. Kokocinski, E. Kulesha, D. London, I. Longden, C. Melsopp, P. Meidl, B. Overduin, A. Parker, G. Proctor, A. Prlic, M. Rae, D. Rios, S. Redmond, M. Schuster, I. Sealy, S. Searle, J. Severin, G. Slater, D. Smedley, J. Smith, A. Stabenau, J. Stalker, S. Trevanion, A. Ureta- Vidal, J. Vogel, S. White, C.Woodwark, and T. J. Hubbard. Ensembl 2006. Nucleic acids research, 34(Database issue), January 2006.
Weigthed total least square regression according to Golub and Van Loan (1980) in SIAM J.Numer.Anal Vol 17 No.6.
tls(formula, D = NULL, T = NULL, precision = .Machine$double.eps)
tls(formula, D = NULL, T = NULL, precision = .Machine$double.eps)
formula |
An object of class formula. |
D |
Diagonal weigth matrix. Default weights are set to 1. |
T |
Diagonal weigth matrix. Default weights are set to 1. |
precision |
Smallest possible numeric value on this machine (default). |
tls
returns a lm object.
Sebastian Duemcke [email protected]
Golub, G.H. and Van Loan, C.F. (1980). An analysis of the total least squares problem. SIAM J. Numer. Anal., 17:883-893.
f = 1.5 # true ratio a = rnorm(5000) b = f*a a = a + rnorm(5000,sd=0.5) b = b + rnorm(5000,sd=0.5) coeff.tls = coef(tls(b ~ a + 0)) coeff.lm1 = coef(lm(b ~ a + 0)) coeff.lm2 = 1/coef(lm(a ~ b + 0)) heatscatter(a,b) abline(0,coeff.lm1,col="red",pch=19,lwd=2) abline(0,coeff.lm2,col="orange",pch=19,lwd=2) abline(0,coeff.tls,col="green",pch=19,lwd=2) abline(0,f,col="grey",pch=19,lwd=2,lty=2) legend("topleft", c("Least-squares regr. (y ~ x + 0)", "Least-squares regr. (x ~ y + 0)", "Total Least-squares regr.", "True ratio"), col=c("red", "orange", "green", "grey"), lty=c(1,1,1,2), lwd=2) results = c(coeff.tls,coeff.lm1,coeff.lm2) names(results) = c("coeff.tls","coeff.lm1","coeff.lm2") print(results)
f = 1.5 # true ratio a = rnorm(5000) b = f*a a = a + rnorm(5000,sd=0.5) b = b + rnorm(5000,sd=0.5) coeff.tls = coef(tls(b ~ a + 0)) coeff.lm1 = coef(lm(b ~ a + 0)) coeff.lm2 = 1/coef(lm(a ~ b + 0)) heatscatter(a,b) abline(0,coeff.lm1,col="red",pch=19,lwd=2) abline(0,coeff.lm2,col="orange",pch=19,lwd=2) abline(0,coeff.tls,col="green",pch=19,lwd=2) abline(0,f,col="grey",pch=19,lwd=2,lty=2) legend("topleft", c("Least-squares regr. (y ~ x + 0)", "Least-squares regr. (x ~ y + 0)", "Total Least-squares regr.", "True ratio"), col=c("red", "orange", "green", "grey"), lty=c(1,1,1,2), lwd=2) results = c(coeff.tls,coeff.lm1,coeff.lm2) names(results) = c("coeff.tls","coeff.lm1","coeff.lm2") print(results)
The phenotype matrix Wt.phenomat
contains information about the experimental design. It is comprised of the filename, the type of RNA fraction measured (T, U or L), the labeling time and the replicate number.
Wt.phenomat
Wt.phenomat
The phenomat is a matrix comprised of the file name, the type of RNA fraction mesasured (T, U or L, fraction
column), the labeling time (time
,timeframe
column) and the replicate number (nr
column). Rows in this matrix represent the individual experiments.
M. Sun, B. Schwalb, D. Schulz, N. Pirkl, L. Lariviere, K. Maier, A. Tresch, P. Cramer. Mutual feedback between mRNA synthesis and degradation buffers transcript levels in a eukaryote. Under review.