knitr::opts_chunk$set(
  collapse = TRUE,
  comment = "#>",
  fig.path = "man/figures/README-",
  out.width = "100%"
)

splice2neo

R-CMD-check Codecov test coverage r badger::badge_devel("TRON-Bioinformatics/splice2neo", "blue") r badger::badge_lifecycle("experimental", "blue") r badger::badge_last_commit("TRON-Bioinformatics/splice2neo")

Documentation: https://tron-bioinformatics.github.io/splice2neo/

Overview

This package provides functions for the analysis of splice junctions and their association with somatic mutations. It integrates the output of several tools that predict splicing effects from mutations or detect expressed splice junctions from RNA-seq data into a standardized splice junction format based on genomic coordinates. Users can filter detected splice junctions against canonical splice junctions and annotate them with affected transcript sequences, CDS, and resulting peptide sequences. Splice2neo supports splice events from alternative 3'/5' splice sites, exons skipping, intron retentions, exitrons, and mutually exclusive exons. Integrating splice2neo functions and detection rules based on splice effect scores and RNA-seq support facilitates the identification of mutation-associated splice junctions that are tumor-specific and can encode neoantigen candidates.

Currently, splice2neo supports the output data from the following tools:

For a more detailed description and a full example workflow see the vignette

Installation

The R package splice2neo is not yet on CRAN or Bioconductor. Therefore, you need to install splice2neo from github.

if (!requireNamespace("remotes", quietly = TRUE)) {
  install.packages("remotes")
}

remotes::install_github("TRON-Bioinformatics/splice2neo")

Example usage

This is a basic example of how to use some functions.

library(splice2neo)

# load human genome reference sequence
requireNamespace("BSgenome.Hsapiens.UCSC.hg19", quietly = TRUE)
bsg <- BSgenome.Hsapiens.UCSC.hg19::BSgenome.Hsapiens.UCSC.hg19

1 Example data

We start with some example splice junctions provided with the package.

junc_df <- dplyr::tibble(
  junc_id = toy_junc_id[c(1, 6, 10)]
)

junc_df

2 Add transcripts

Next, we find the transcripts which are in the same genomic region as the splice junction and that may be affected by the junction.

junc_df %>% 
  add_tx(toy_transcripts) %>% 
  head()

3 Modify transcripts with junctions

We modify the canonical transcripts by introducing the splice junctions. Then we add the transcript sequence in a fixed-sized window around the junction positions, the context sequence.

toy_junc_df %>% 
  head()


toy_junc_df %>% 
  add_context_seq(transcripts = toy_transcripts, size = 400, bsg = bsg) %>% 
  head()

4 Annotate peptide sequence

Here, we use the splice junctions to modify the coding sequences (CDS) of the reference transcripts. The resulting CDS sequences are translated into protein sequence and further annotated with the peptide around the junction, the relative position of the splice junction in the peptide, and the location of the junction in an open reading frame (ORF).

toy_junc_df %>% 

  # add peptide sequence
  add_peptide(cds=toy_cds, flanking_size = 13, bsg = bsg) %>% 

  # select subset of columns
  dplyr::select(junc_id, peptide_context, junc_in_orf, cds_description) %>% 
  head()

Issues

Please report issues here: https://github.com/TRON-Bioinformatics/splice2neo/issues

Citation

Lang, Franziska, Patrick Sorn, Martin Suchan, Alina Henrich, Christian Albrecht, Nina Koehl, Aline Beicht, et al. 2023. “Prediction of Tumor-Specific Splicing from Somatic Mutations as a Source of Neoantigen Candidates.” bioRxiv. https://doi.org/10.1101/2023.06.27.546494.



TRON-Bioinformatics/splice2neo documentation built on July 1, 2024, 7:57 p.m.