R/pStrataTEI.R
pStrataTEI.Rd
This function computes the partial transcriptome evolutionary index (TEI) values combined for each strata.
In detail, each gene gets a TEI contribution profile.
$$TEI_is = f_is * ps_i$$
where TEI_is is the partial TEI value of gene i, \(f_is = e_is / \sum e_is\) and \(ps_i\) is the phylostratum of gene i.
\(TEI_is\) values are combined per \(ps\).
pStrataTEI(
ExpressionSet,
Phylostratum = NULL,
by = NULL,
groups = NULL,
split = 1e+05,
showprogress = TRUE,
threads = 1
)
expression object with rownames as GeneID (dgCMatrix) or standard PhyloExpressionSet object.
a named vector representing phylostratum per GeneID with names as GeneID (not used if Expression is PhyloExpressionSet).
specify min/max transformation by row (stages, cells, celltypes) or by column (phylostratum)
specify stages or cells to be grouped into celltypes by a named list
specify number of columns to split
boolean if progressbar should be shown
specify number of threads
a numeric matrix storing the summed partial TEI values for each strata.
The partial TEI values combined per strata give an overall
impression of the contribution of each
strata to the global TEI
pattern.
Domazet-Loso T. and Tautz D. (2010). A phylogenetically based transcriptome age index mirrors ontogenetic divergence patterns. Nature (468): 815-818.
Quint M et al. (2012). A transcriptomic hourglass in plant embryogenesis. Nature (490): 98-101.
Drost HG et al. (2015) Mol Biol Evol. 32 (5): 1221-1231 doi:10.1093/molbev/msv012
## get Seurat object
celegans<-readRDS(file=system.file("extdata",
"celegans.embryo.SeuratData.rds", package="scTEI")
)
## load Caenorhabditis elegans gene age estimation
celegans_ps<-readr::read_tsv(
file=system.file("extdata",
"Sun2021_Orthomap.tsv", package="scTEI")
)
#> Rows: 20040 Columns: 2
#> ── Column specification ────────────────────────────────────────────────────────
#> Delimiter: "\t"
#> chr (1): GeneID
#> dbl (1): Phylostratum
#>
#> ℹ Use `spec()` to retrieve the full column specification for this data.
#> ℹ Specify the column types or set `show_col_types = FALSE` to quiet this message.
## define Phylostratum
ps_vec<-setNames(
as.numeric(celegans_ps$Phylostratum),
celegans_ps$GeneID
)
## get partial TEI strata values
Seurat::Idents(celegans)<-"embryo.time.bin"
pSt<-pStrataTEI(
ExpressionSet=Seurat::GetAssayData(celegans, assay="RNA", layer="counts"),
Phylostratum=ps_vec
)
## get partial TEI strata values per cell group
Seurat::Idents(celegans)<-"embryo.time.bin"
cell_groups<-Ident2cellList(Seurat::Idents(celegans))
pSt<-pStrataTEI(
ExpressionSet=Seurat::GetAssayData(celegans, assay="RNA", layer="counts"),
Phylostratum=ps_vec,
groups=cell_groups
)
p1<-ComplexHeatmap::Heatmap(
pSt,
name="pSt",
column_title="partial TEI values - cell groups",
row_title="More Recent <<< More Ancient",
cluster_rows=FALSE,
cluster_columns=FALSE,
col=viridis::viridis(3)
)
p1
## get relative expression over stages per cell group
Seurat::Idents(celegans)<-"embryo.time.bin"
cell_groups<-Ident2cellList(Seurat::Idents(celegans))
pSt<-pStrataTEI(
ExpressionSet=Seurat::GetAssayData(celegans, assay="RNA", layer="counts"),
Phylostratum=ps_vec,
groups=cell_groups,
by="row"
)
p2<-ComplexHeatmap::Heatmap(
pSt,
name="pSt",
column_title="partial TEI values - cell groups - by row",
row_title="More Recent <<< More Ancient",
cluster_rows=FALSE,
cluster_columns=FALSE,
col=viridis::viridis(3)
)
p2
## get relative expression over phylostrata per cell group
Seurat::Idents(celegans)<-"embryo.time.bin"
cell_groups<-Ident2cellList(Seurat::Idents(celegans))
pSt<-pStrataTEI(
ExpressionSet=Seurat::GetAssayData(celegans, assay="RNA", layer="counts"),
Phylostratum=ps_vec,
groups=cell_groups,
by="column"
)
p3<-ComplexHeatmap::Heatmap(
pSt,
name="pSt",
column_title="partial TEI values - cell groups - by col",
row_title="More Recent <<< More Ancient",
cluster_rows=FALSE,
cluster_columns=FALSE,
col=viridis::viridis(3)
)
p3