Differential_genes: Calculate differential genes between cells/clusters.

View source: R/GroupsDiffGenes.R

Differential_genesR Documentation

Calculate differential genes between cells/clusters.

Description

Provide differential genes between given groups.

Usage

Differential_genes(
  data_list = list(),
  min_Sample = 5,
  min_Gene = 1500,
  DDLK_Clusters,
  data_id,
  Genesets,
  data_type = c("Normalised", "Raw"),
  DifferentiateBy = "Clusters",
  p_val = 0.05,
  lfc = 0,
  up_gene_number = 10
)

Arguments

data_list

List of expression matricies

min_Sample

gene filter, filter out genes which are not expressed in at least min_Sample cells

min_Gene

cell filter, filter out those cells which do not express at least min_Gene genes

DDLK_Clusters

Output of DDLK_Clust.R method

data_id

List of names/ids of expression matrix

Genesets

list of genesets/pathways used to calculate PathwayEnrichmentScore.

data_type

list of expression data passed in data. Valid inputs are either raw or normalised.

DifferentiateBy

Any column name from DDLK_Clusters$PathwayDDLK_clust, default is "Clusters"

p_val

Threshold p value, default is 0.05

lfc

Threshold log fold change value, default is 0

up_gene_number

select number of upregulated genes in each group

Value

Diff_mat list of three objects 1. DifferentialMatrix = Data matrix with top selected up_gene_number in each group 2. Diffup_genes = list of differential gene in each group with selected up_gene_number 3. annotations = Cell wise annotation of DifferentialMatrix

Examples

data1 = unCTC::Poonia_et_al._TPMData
data2 = unCTC::Ding_et_al._WBC1_TPMData
Data_list = list(data1,data2)
Data_Id = list("data1","data2")
Genesets = unCTC::c2.all.v7.2.symbols
Pathway_score = PathwayEnrichmentScore(data_list=Data_list,
                                        data_id= Data_Id,
                                         min_Sample = 5,
                                         min_Gene = 1500,
                                        Genesets=Genesets,
                                        min.size=70,
                                        max.size=100)

DDLK_Clusters = DDLK_Clust(PathwayScore = Pathway_score$Pathway_score,
                           PathwayMetaData=Pathway_score$Pathway_metadata,
                            n=3,
                            out.dir = paste0(getwd(),"/unCTC"))

Output = Differential_genes(data_list=Data_list,
                              min_Sample = 5,
                              min_Gene = 1500,
                              DDLK_Clusters=DDLK_Clusters,
                              data_id = Data_Id,
                              Genesets=Genesets,
                              data_type="Normalised",
                              DifferentiateBy = "Clusters")

SaritaPoonia/unCTC documentation built on Nov. 8, 2022, 12:07 p.m.