calc_diff_abund_deseq2: Differential abundance with DESeq2

View source: R/calculations--differential_abundance.R

calc_diff_abund_deseq2R Documentation

Differential abundance with DESeq2

Description

EXPERIMENTAL: This function is still being tested and developed; use with caution. Uses the DESeq2-package package to conduct differential abundance analysis of count data. Counts can be of OTUs/ASVs or taxa. The plotting function heat_tree_matrix is useful for visualizing these results. See details section below for considerations on preparing data for this analysis.

Usage

calc_diff_abund_deseq2(
  obj,
  data,
  cols,
  groups,
  other_cols = FALSE,
  lfc_shrinkage = c("none", "normal", "ashr"),
  ...
)

Arguments

obj

A taxmap object

data

The name of a table in obj that contains data for each sample in columns.

cols

The names/indexes of columns in data to use. By default, all numeric columns are used. Takes one of the following inputs:

TRUE/FALSE:

All/No columns will used.

Character vector:

The names of columns to use

Numeric vector:

The indexes of columns to use

Vector of TRUE/FALSE of length equal to the number of columns:

Use the columns corresponding to TRUE values.

groups

A vector defining how samples are grouped into "treatments". Must be the same order and length as cols.

other_cols

If TRUE, preserve all columns not in cols in the output. If FALSE, dont keep other columns. If a column names or indexes are supplied, only preserve those columns.

lfc_shrinkage

What technique to use to adjust the log fold change results for low counts. Useful for ranking and visualizing log fold changes. Must be one of the following:

'none'

No log fold change adjustments.

'normal'

The original DESeq2 shrinkage estimator

'ashr'

Adaptive shrinkage estimator from the ashr package, using a fitted mixture of normals prior.

...

Passed to results if the lfc_shrinkage option is "none" and to lfcShrink otherwise.

Details

Data should be raw read counts, not rarefied, converted to proportions, or modified with any other technique designed to correct for sample size since DESeq2-package is designed to be used with count data and takes into account unequal sample size when determining differential abundance. Warnings will be given if the data is not integers or all sample sizes are equal.

Value

A tibble with at least the taxon ID of the thing tested, the groups compared, and the DESeq2 results. The log2FoldChange values will be positive if treatment_1 is more abundant and treatment_2.

See Also

Other calculations: calc_group_mean(), calc_group_median(), calc_group_rsd(), calc_group_stat(), calc_n_samples(), calc_obs_props(), calc_prop_samples(), calc_taxon_abund(), compare_groups(), counts_to_presence(), rarefy_obs(), zero_low_counts()

Examples

## Not run: 
# Parse data for plotting
x = parse_tax_data(hmp_otus, class_cols = "lineage", class_sep = ";",
                   class_key = c(tax_rank = "taxon_rank", tax_name = "taxon_name"),
                   class_regex = "^(.+)__(.+)$")

# Get per-taxon counts
x$data$tax_table <- calc_taxon_abund(x, data = "tax_data", cols = hmp_samples$sample_id)

# Calculate difference between groups
x$data$diff_table <- calc_diff_abund_deseq2(x, data = "tax_table",
                                    cols = hmp_samples$sample_id,
                                    groups = hmp_samples$body_site)
                                    
# Plot results (might take a few minutes)
heat_tree_matrix(x,
                 data = "diff_table",
                 node_size = n_obs,
                 node_label = taxon_names,
                 node_color = ifelse(is.na(padj) | padj > 0.05, 0, log2FoldChange),
                 node_color_range = diverging_palette(),
                 node_color_trans = "linear",
                 node_color_interval = c(-3, 3),
                 edge_color_interval = c(-3, 3),
                 node_size_axis_label = "Number of OTUs",
                 node_color_axis_label = "Log2 fold change")


## End(Not run)


metacoder documentation built on April 4, 2023, 9:08 a.m.