Identifies transcripts/genes that are consistently expressed above a threshold across samples. This function adds a logical column `.abundant` to indicate which features pass the filtering criteria.

identify_abundant(
  .data,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  factor_of_interest = NULL,
  design = NULL,
  minimum_counts = 10,
  minimum_proportion = 0.7
)

# S4 method for class 'spec_tbl_df'
identify_abundant(
  .data,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  factor_of_interest = NULL,
  design = NULL,
  minimum_counts = 10,
  minimum_proportion = 0.7
)

# S4 method for class 'tbl_df'
identify_abundant(
  .data,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  factor_of_interest = NULL,
  design = NULL,
  minimum_counts = 10,
  minimum_proportion = 0.7
)

# S4 method for class 'tidybulk'
identify_abundant(
  .data,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  factor_of_interest = NULL,
  design = NULL,
  minimum_counts = 10,
  minimum_proportion = 0.7
)

# S4 method for class 'SummarizedExperiment'
identify_abundant(
  .data,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  factor_of_interest = NULL,
  design = NULL,
  minimum_counts = 10,
  minimum_proportion = 0.7
)

# S4 method for class 'RangedSummarizedExperiment'
identify_abundant(
  .data,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  factor_of_interest = NULL,
  design = NULL,
  minimum_counts = 10,
  minimum_proportion = 0.7
)

Arguments

.data

A `tbl` or `SummarizedExperiment` object containing transcript/gene abundance data

.sample

The name of the sample column

.transcript

The name of the transcript/gene column

.abundance

The name of the transcript/gene abundance column

factor_of_interest

The name of the column containing groups/conditions for filtering. Used by edgeR's filterByExpr to define sample groups.

design

A design matrix for more complex experimental designs. If provided, this is passed to filterByExpr instead of factor_of_interest.

minimum_counts

A positive number specifying the minimum counts per million (CPM) threshold for a transcript to be considered abundant (default = 10)

minimum_proportion

A number between 0 and 1 specifying the minimum proportion of samples that must exceed the minimum_counts threshold (default = 0.7)

Value

Returns the input object with an additional logical column `.abundant` indicating which features passed the abundance threshold criteria.

A consistent object (to the input) with additional columns for the statistics from the hypothesis test (e.g., log fold change, p-value and false discovery rate).

A consistent object (to the input) with additional columns for the statistics from the hypothesis test (e.g., log fold change, p-value and false discovery rate).

A consistent object (to the input) with additional columns for the statistics from the hypothesis test (e.g., log fold change, p-value and false discovery rate).

A `SummarizedExperiment` object

A `SummarizedExperiment` object

Details

`r lifecycle::badge("maturing")`

This function uses edgeR's filterByExpr() function to identify consistently expressed features. A feature is considered abundant if it has CPM > minimum_counts in at least minimum_proportion of samples in at least one experimental group (defined by factor_of_interest or design).

References

McCarthy, D. J., Chen, Y., & Smyth, G. K. (2012). Differential expression analysis of multifactor RNA-Seq experiments with respect to biological variation. Nucleic Acids Research, 40(10), 4288-4297. DOI: 10.1093/bioinformatics/btp616

Examples

# Basic usage
se_mini |> identify_abundant()
#> No group or design set. Assuming all samples belong to one group.
#> class: SummarizedExperiment 
#> dim: 527 5 
#> metadata(0):
#> assays(1): count
#> rownames(527): ABCB4 ABCB9 ... ZNF324 ZNF442
#> rowData names(2): entrez .abundant
#> colnames(5): SRR1740034 SRR1740035 SRR1740043 SRR1740058 SRR1740067
#> colData names(5): Cell.type time condition days dead

# With custom thresholds
se_mini |> identify_abundant(
  minimum_counts = 5,
  minimum_proportion = 0.5
)
#> No group or design set. Assuming all samples belong to one group.
#> class: SummarizedExperiment 
#> dim: 527 5 
#> metadata(0):
#> assays(1): count
#> rownames(527): ABCB4 ABCB9 ... ZNF324 ZNF442
#> rowData names(2): entrez .abundant
#> colnames(5): SRR1740034 SRR1740035 SRR1740043 SRR1740058 SRR1740067
#> colData names(5): Cell.type time condition days dead

# Using a factor of interest
se_mini |> identify_abundant(factor_of_interest = condition)
#> class: SummarizedExperiment 
#> dim: 527 5 
#> metadata(0):
#> assays(1): count
#> rownames(527): ABCB4 ABCB9 ... ZNF324 ZNF442
#> rowData names(2): entrez .abundant
#> colnames(5): SRR1740034 SRR1740035 SRR1740043 SRR1740058 SRR1740067
#> colData names(5): Cell.type time condition days dead