ggmRSCU {ggmRSCU}R Documentation

Visualization of RSCU (Relative Synonymous Codon Usage) with ggplot2

Description

This function generates a customized visualization of Relative Synonymous Codon Usage (RSCU) across distinct biological categories. It represents codon usage through color-coded blocks using species-specific shapes, and optionally incorporates a secondary plot displaying codon labels. Also applicable for extended visualizations of similar data structures (e.g., energy consumption data).

Usage

ggmRSCU(
  data,
  x,
  y,
  fill,
  rev = FALSE,
  theme = "theme1",
  shape = FALSE,
  size = 1,
  width = 0.1,
  sub_axis = "fill",
  merge = FALSE,
  border_color = NA,
  blocks_colors = NULL
)

Arguments

data

A data frame containing the input data.

x

Primary and secondary x-axis variables: a character vector (e.g., c("Group", "Subgroup")) or a single variable ("Group").

y

Primary and secondary y-axis variables: a character vector (e.g., c("ValueVar", "Subgroup")) or a single variable ("ValueVar").

fill

Variable for color filling.

rev

(Optional) logical value. If TRUE, the color order will be reversed. The default is FALSE.

theme

(Optional) A string specifying the color theme for the blocks. Available themes are: "theme1"-"theme8".

shape

(Optional) logical value. If TRUE, a vector of shape values for species in the plot is applied (default: FALSE).

size

(Optional) Size of shape markers (default: 1).

width

(Optional) Thickness of bar plot borders. (e.g., 0.1 for thin, 0.5 for thick). (Default: 0.1).

sub_axis

(Optional) Annotation type for secondary plot ("ysub" or "fill", default: "fill").

merge

(Optional) logical value. If TRUE, the main plot will be merged with the secondary annotation plot. Default is FALSE.

border_color

(Optional) Color for bar borders (default: NA)

blocks_colors

(Optional) A custom color vector to fill the plot. If "NULL", the default color scheme will be used.

Details

Key features: - Hierarchical grouping with automatic position calculation - Flexible color management with 8 predefined themes or custom palettes - Dual-axis annotation system for biological information - Interactive subplot arrangements with smart label positioning

Value

A ggplot object or patchwork composite plot when merge=TRUE

Examples

# Example 1: Basic Usage

library(ggmRSCU)

ggmRSCU(
  data = rscu_data,
  x = c(AA, Species),
  y = c(RSCU, Codon),
  fill = Fill,
  shape = TRUE,
  merge = TRUE
)


# Example 2: Visualize RSCU from data loaded using `read_codonw` function

# Load data from `read_codonw` output
example_dir <- system.file("extdata", "codonw", package = "ggmRSCU")
data <- read_codonw(example_dir)

col <- c("#FF6F61", "#6B5B95", "#88B04B",
         "#F7CAC9", "#92A8D1", "#9b59b6")

ggmRSCU(
  data = data,
  x = c(AA, Species),
  y = c(RSCU, Codon),
  fill = Fill,
  blocks_colors = col,
  shape = TRUE,
  merge = TRUE,
  sub_axis = "fill"
)



# Example 3: Visualize RSCU from data loaded using `read_mega` function

# Load data from `read_mega` output
example_dir <- system.file("extdata", "mega", package = "ggmRSCU")
df <- read_mega(example_dir)

ggmRSCU(
  data = df,
  x = c(AA, Species),
  y = c(RSCU, Codon),
  fill = Fill
)



# Example 4: Additional Dataset (scoring_data)

col <- c("#4E79A7", "#A0CBE8", "#F28E2B",
         "#FFBE7D", "#59A14F", "#9b59b6")

ggmRSCU(
  data = scoring_data,
  x = c(Sample, Replicate),
  y = c(Score, CellType),
  fill = CellType,
  blocks_colors = col,
  shape = TRUE,
  sub_axis = "Fill"
)


# Example 5: Additional Dataset (hydrogen_data)

library(tidyr)

df_long <- hydrogen_data %>%
  pivot_longer(cols = -c(Year, type),
  names_to = "category",
  values_to = "value")

col <- c("#4E79A7", "#A0CBE8", "#F28E2B",
         "#FFBE7D", "#59A14F", "#9b59b6")

ggmRSCU(
  data = df_long,
  x = c(Year, category),
  y = c(value, type),
  fill = type,
  shape = TRUE,
  blocks_colors = col
)


# Example 6: Additional Dataset (energy_data)

long_data <- energy_data %>%
  pivot_longer(cols = 3:7,
  names_to = "Source",
  values_to = "Production")

col <- c("#4E79A7", "#A0CBE8", "#F28E2B",
         "#FFBE7D", "#59A14F", "#9b59b6")

ggmRSCU(
  data = long_data,
  x = c(Scenario, Year),
  y = c(Production, Source),
  fill = Source,
  blocks_colors = col,
  shape = TRUE
)


# Example 7: Additional Dataset (financial_data)

ggmRSCU(
  data = financial_data,
  x = c(Month, Category),
  y = c(Amount),
  fill = Subcategory
)


# Example 8: Additional Dataset (genomic_data)

ggmRSCU(
  data = genomic_data,
  x = c(Species, haplotypes),
  y = c(Mb),
  fill = Region,
  shape = TRUE
)



[Package ggmRSCU version 0.1.0 Index]