Title: | Differential Analysis of Rhythmic Transcriptome Data |
---|---|
Description: | A flexible approach, inspired by cosinor regression, for differential analysis of rhythmic transcriptome data. See Singer and Hughey (2018) <doi:10.1177/0748730418813785>. |
Authors: | Jake Hughey [aut, cre], Jordan Singer [ctb] |
Maintainer: | Jake Hughey <[email protected]> |
License: | GPL-2 |
Version: | 1.0.2 |
Built: | 2025-01-19 03:46:12 UTC |
Source: | https://github.com/hugheylab/limorhyde |
Generate basis matrix for cosinor regression.
getCosinorBasis(x, period, intercept)
getCosinorBasis(x, period, intercept)
x |
Values of the predictor variable. |
period |
Period for the predictor variable. |
intercept |
If |
A matrix with a row for each value of x
and a column for each
component of the decomposition.
b = getCosinorBasis(seq(0, 20, 4), period = 24, intercept = FALSE)
b = getCosinorBasis(seq(0, 20, 4), period = 24, intercept = FALSE)
Generate basis matrix for a periodic B-spline using pbs::pbs()
.
getSplineBasis(x, period, nKnots, intercept)
getSplineBasis(x, period, nKnots, intercept)
x |
Values of the predictor variable. |
period |
Period for the predictor variable. |
nKnots |
Number of internal knots. |
intercept |
If |
A matrix with a row for each value of x
and a column for each
component of the decomposition.
b = getSplineBasis(seq(0, 20, 4), period = 24, nKnots = 3, intercept = FALSE)
b = getSplineBasis(seq(0, 20, 4), period = 24, nKnots = 3, intercept = FALSE)
Decompose a periodic time variable into multiple components based on either the first harmonic of a Fourier series or on a periodic smoothing spline.
limorhyde( time, colnamePrefix = NULL, period = 24, sinusoid = TRUE, nKnots = 3, intercept = FALSE )
limorhyde( time, colnamePrefix = NULL, period = 24, sinusoid = TRUE, nKnots = 3, intercept = FALSE )
time |
Numeric vector of times, e.g., at which samples were acquired. |
colnamePrefix |
Character string with which to prefix the column names of the basis. |
period |
Number corresponding to the period to use for the
decomposition (in the same units as |
sinusoid |
If |
nKnots |
Number of internal knots for the periodic spline. Only used if
|
intercept |
If |
A matrix with a row for each sample and a column for each component of the time decomposition.
# create an example data frame nSamples = 12 d = data.frame( sample = paste0('sample_', 1:nSamples), genotype = factor(rep(c('WT', 'KO'), each = nSamples / 2), levels = c('WT', 'KO')), zt = rep(seq(0, 24 - 24 / nSamples * 2, 24 / nSamples * 2), times = 2), stringsAsFactors = FALSE) # call limorhyde limo = limorhyde(d$zt, 'zt_') d = cbind(d, limo) # create a design matrix that could be used with methods such as limma design = model.matrix(~ genotype * (zt_cos + zt_sin), data = d)
# create an example data frame nSamples = 12 d = data.frame( sample = paste0('sample_', 1:nSamples), genotype = factor(rep(c('WT', 'KO'), each = nSamples / 2), levels = c('WT', 'KO')), zt = rep(seq(0, 24 - 24 / nSamples * 2, 24 / nSamples * 2), times = 2), stringsAsFactors = FALSE) # call limorhyde limo = limorhyde(d$zt, 'zt_') d = cbind(d, limo) # create a design matrix that could be used with methods such as limma design = model.matrix(~ genotype * (zt_cos + zt_sin), data = d)