[BioC] multiple groups time course RNA Seq LIMMA

Gordon K Smyth smyth at wehi.EDU.AU
Fri Mar 14 01:18:00 CET 2014


Dear Michela Riba,

The first thing that puzzles me that the design matrix you say you have 
used would give no residual degrees of freedom for your data.  In this 
situation, it would impossible for you to do any testing for DE genes, and 
limma would have given you a warning about this.  Given that you have only 
5 time points, the largest df you can use for the spline curves is 3, and 
you might consider using just 2 df.

Anyway, I will go back to the original question.  You haven't said 
anything about replicates, but I assume that you have unreplicated time 
courses.  There two main questions you can ask.  You can find genes that 
are DE within each time course separately, or you can find genes that 
differ between time courses.

You have 6 groups and 5 times points.  Here is a simulated version with 
1000 genes:

   Group <- gl(6,5,30)
   Time <- as.numeric(gl(5,1,30))
   X <- ns(Time,df=2)
   y <- matrix(rnorm(1000*30),1000,30)

FIRST, how to find genes that are DE within each time course:

   design <- model.matrix(~Group+Group:X)
   colnames(design) <- make.names(colnames(design))
   fit <- lmFit(y,design)
   fit <- eBayes(fit)

This gives:

  > colnames(fit)
   [1] "X.Intercept." "Group2"       "Group3"       "Group4"
   [5] "Group5"       "Group6"       "Group1.X1"    "Group2.X1"
   [9] "Group3.X1"    "Group4.X1"    "Group5.X1"    "Group6.X1"
  [13] "Group1.X2"    "Group2.X2"    "Group3.X2"    "Group4.X2"
  [17] "Group5.X2"    "Group6.X2"

To find genes changing over time within the first group (control):

   topTable(fit,coef=c("Group1.X1","Group1.X2"))

To find genes changing over time in the second group:

   topTable(fit,coef=c("Group2.X1","Group2.X2"))

And so on.  You get the idea?

SECOND, how to find genes that have different responses between the 
different time courses.

To compare Group2 with Group1:

   cont.mat <- makeContrasts(
      Group2.X1-Group1.X1,
      Group2.X2-Group1.X2,
      levels=design)
   fit2 <- contrasts.fit(fit,cont.mat)
   fit2 <- eBayes(fit2)
   topTable(fit2)

To compare any group with any other, you only need to change the group 
numbers in the makeContrasts command.

Best wishes
Gordon


On Wed, 12 Mar 2014, Riba Michela wrote:
...
> the design is derived from the limma vignette
>
> X <- ns(limmaTargets$Time,df=4)
> design<-model.matrix(~Group*X) #here Group refers to the different 6 experimental groups)
>
> My point regards the possibility and way to extract genes which differ 
> in the different experimental groups (that unfortunately yesterday I 
> mentioned as treatments) not compared with the control, or better which 
> could be the best way to solve this, for example reasoning simply on 
> venn diagram of the single results as was for example done using the 
> maSigPro package for time series analysis in gene expression.
>
> The point regards exactly which could be the most convenient way to 
> design the analysis
> -compare just 2 time series (control and experimental group1) just with 
> control and one experimental group
> -make a design matrix to compare different experimental groups with 
> control o among each other
...


> Date: Tue, 14 Jan 2014 11:53:11 +0100
> From: Riba Michela <riba.michela at hsr.it>
> To: "bioconductor at r-project.org" <bioconductor at r-project.org>
> Subject: [BioC] multiple groups time course RNA Seq LIMMA
>
>
> Hi,

> I'm approaching a RNA-seq experiment concerning the analysis of a time 
course of 5 time points in 6 experimental groups (including Control
> group).

> As an example:
>
> FileName Group Time
> a Control 6hr
> b Control 24h
> ...
> e ExpG1 6hr
> f ExpG1 24hr
> ...
> l ExpG2 6hr
> m ExpG2 24hr
> ...
> (ExpG1, ExpG2 are experimental groups)
>
> I'm going to use LIMMA for extraction of time changing genes in the 
> single experimental groups compared to Control group.
>
> I'd like to see how to extract this result in the topTable (i.e. which 
> coefficients select) for each single comparison of the experimental 
> group towards Control) since from the provided example in the LIMMA 
> manual (pag. 49) such topTable is referred to a design concerning one 
> single experimental group towards Control in a time course instead of 
> multiple experimental groups.
>
> Is there in addition the possibility to design a contrast matrix in such 
> situations or is it better to consider topTables using various
> coefficients blocks?
>
>
> Thanks a lot for your answer
> Michela Riba

______________________________________________________________________
The information in this email is confidential and intend...{{dropped:4}}



More information about the Bioconductor mailing list