
Generate Beta Diversity Ordination Plot for Cross-Sectional Data
Source:R/generate_beta_ordination_single.R
generate_beta_ordination_single.RdCreates PCoA ordination plots for single time point or cross-sectional microbiome data with optional group comparisons.
Usage
generate_beta_ordination_single(
data.obj,
time.var = NULL,
t.level = NULL,
group.var = NULL,
adj.vars = NULL,
strata.var = NULL,
dist.obj = NULL,
dist.name = c("BC", "Jaccard"),
pc.obj = NULL,
base.size = 16,
point.size = 4,
theme.choice = "bw",
custom.theme = NULL,
palette = NULL,
pdf = TRUE,
file.ann = NULL,
pdf.wid = 11,
pdf.hei = 8.5,
...
)Arguments
- data.obj
A MicrobiomeStat data object, which is a list containing at minimum the following components:
feature.tab: A matrix of feature abundances (taxa/genes as rows, samples as columns)meta.dat: A data frame of sample metadata (samples as rows)
Optional components include:
feature.ann: A matrix/data frame of feature annotations (e.g., taxonomy)tree: A phylogenetic tree object (class "phylo")feature.agg.list: Pre-aggregated feature tables by taxonomy
Data objects can be created using converters like
mStat_convert_phyloseq_to_data_objor importers likemStat_import_qiime2_as_data_obj.- time.var
Character string specifying the column name in meta.dat containing the time variable. Required for longitudinal and paired analyses. Should be a factor or character with meaningful time point labels.
- t.level
Character string specifying the time level to subset to. If NULL, uses all data.
- group.var
Character string specifying the column name in meta.dat containing the grouping variable (e.g., treatment, condition, phenotype). Used for between-group comparisons.
- adj.vars
Character vector specifying column names in meta.dat to be used as covariates for adjustment in statistical models. These variables will be included as fixed effects.
- strata.var
Character string specifying the column name in meta.dat for stratification. When provided, analyses and visualizations will be performed separately within each stratum (e.g., by site, batch, or sex).
- dist.obj
A list of pre-calculated distance matrices. If NULL and distances are needed, they will be calculated automatically. List names should match dist.name (e.g., "BC" for Bray-Curtis). See
mStat_calculate_beta_diversity.- dist.name
Character vector specifying which distance metrics to use. Options depend on available methods:
"BC": Bray-Curtis dissimilarity
"Jaccard": Jaccard distance
"UniFrac": Unweighted UniFrac (requires tree)
"GUniFrac": Generalized UniFrac (requires tree)
"WUniFrac": Weighted UniFrac (requires tree)
"JS": Jensen-Shannon divergence
- pc.obj
A list containing dimension reduction results from
mStat_calculate_PC. If NULL, PCoA is performed automatically.- base.size
Numeric value specifying the base font size for plot text elements. Default is typically 16.
- point.size
Numeric. Size of points in the scatter plot. Default is 4.
- theme.choice
Character string specifying the ggplot2 theme to use. Options include:
"bw": Black and white theme (theme_bw)
"classic": Classic theme (theme_classic)
"minimal": Minimal theme (theme_minimal)
"prism": GraphPad Prism-like theme
"nature": Nature journal style
"light": Light theme (theme_light)
Can also use a custom ggplot2 theme object via custom.theme.
- custom.theme
A custom ggplot2 theme object to override theme.choice. Should be created using ggplot2::theme() or a complete theme function.
- palette
Character vector of colors or a named palette for the plot. If NULL, uses default MicrobiomeStat color scheme. Can be:
A vector of color codes (e.g., c("#E41A1C", "#377EB8"))
A palette name recognized by the plotting function
Logical. If TRUE, saves the plot(s) to PDF file(s) in the current working directory. Default is TRUE.
- file.ann
Character string for additional annotation to append to output filenames. Useful for distinguishing multiple outputs.
- pdf.wid
Numeric value specifying the width of PDF output in inches. Default is typically 11.
- pdf.hei
Numeric value specifying the height of PDF output in inches. Default is typically 8.5.
- ...
Additional arguments passed to underlying functions.
Author
Caffery Yang [email protected]
Examples
if (FALSE) { # \dontrun{
library(aplot)
data(peerj32.obj)
dist.obj <- mStat_calculate_beta_diversity(peerj32.obj, dist.name = c('BC', 'Jaccard'))
pc.obj <- mStat_calculate_PC(dist.obj, method = c('mds'), k = 2, dist.name = c('BC','Jaccard'))
generate_beta_ordination_single(
data.obj = peerj32.obj,
dist.obj = NULL,
pc.obj = NULL,
time.var = "time",
t.level = "2",
group.var = "group",
strata.var = "sex",
adj.vars = "sex",
dist.name = c("BC"),
base.size = 20,
point.size = 4,
theme.choice = "bw",
custom.theme = NULL,
palette = NULL,
pdf = TRUE,
file.ann = NULL,
pdf.wid = 11,
pdf.hei = 8.5
)
data("subset_T2D.obj")
dist.obj <- mStat_calculate_beta_diversity(subset_T2D.obj, dist.name = c('BC', 'Jaccard'))
pc.obj <- mStat_calculate_PC(dist.obj, method = c('mds'), k = 2, dist.name = c('BC','Jaccard'))
generate_beta_ordination_single(
data.obj = subset_T2D.obj,
dist.obj = dist.obj,
pc.obj = pc.obj,
time.var = "visit_number_num",
t.level = NULL,
group.var = "subject_race",
strata.var = "subject_gender",
adj.vars = "sample_body_site",
dist.name = c("BC", 'Jaccard'),
base.size = 20,
point.size = 4,
theme.choice = "bw",
custom.theme = NULL,
palette = NULL,
pdf = TRUE,
file.ann = NULL,
pdf.wid = 11,
pdf.hei = 8.5
)
} # }