Title: | Reversion Mutation Identifier for Sequencing Data |
---|---|
Description: | A tool for detecting reversions for a given pathogenic mutation from next-generation DNA sequencing data. It analyses reads aligned to the locus of the pathogenic mutation and reports reversion events where secondary mutations have restored or undone the deleterious effect of the original pathogenic mutation, e.g., secondary indels complement to a frameshift pathogenic mutation converting the orignal frameshift mutation into inframe mutaions, deletions or SNVs that replaced the original pathogenic mutation restoring the open reading frame, SNVs changing the stop codon caused by the original nonsense SNV into an amino acid, etc. |
Authors: | Hui Xiao [aut, cre], Adam Mills [aut], John Alexander [ctb], Stephen Pettitt [aut], Syed Haider [aut] |
Maintainer: | Hui Xiao <[email protected]> |
License: | GPL-2 |
Version: | 0.0.1 |
Built: | 2024-12-18 06:59:22 UTC |
Source: | CRAN |
Function for detecting reversions for a given pathogenic mutation from reads alignment of NGS genomic data
getReversions( bam.file, genome.version = "BSgenome.Hsapiens.UCSC.hg38", chromosome, pathog.mut.start, pathog.mut.type = "SNV", snv.reference.allele = NULL, snv.alternative.allele = NULL, deletion.sequence = NULL, deletion.length = NULL, insertion.sequence = NULL, flanking.window = 100, minus.strand = FALSE, check.wildtype.reads = FALSE )
getReversions( bam.file, genome.version = "BSgenome.Hsapiens.UCSC.hg38", chromosome, pathog.mut.start, pathog.mut.type = "SNV", snv.reference.allele = NULL, snv.alternative.allele = NULL, deletion.sequence = NULL, deletion.length = NULL, insertion.sequence = NULL, flanking.window = 100, minus.strand = FALSE, check.wildtype.reads = FALSE )
bam.file |
Character. The input bam file containing aligned reads. |
genome.version |
Character. Genome version of alignments in bam file and the position of pathogenic mutation. The name of genome version is specified by the available BSgenome data packages: BSgenome::available.genomes(). Default is "BSgenome.Hsapiens.UCSC.hg38". |
chromosome |
Character. Name of chromosome where pathogenic mutation is located, e.g., "chr17" or "17", "chrX" or "X". The chromosome name should be concordant with the chromosome identifiers used in bam.file and genome.version |
pathog.mut.start |
Integer. Genomic coordinate of the start position of pathogenic mutation. |
pathog.mut.type |
Character. Type of pathogenic mutation: "SNV", "DEL" or "INS". |
snv.reference.allele |
Character. Reference allele of pathogenic mutation if it is a single nucleotide variant (SNV). Default is NULL. |
snv.alternative.allele |
Character. Alternative allele of pathogenic mutation if it is a SNV. Default is NULL. |
deletion.sequence |
Character. Deleted nucleotides of pathogenic mutation if it is a deletion (DEL). Default is NULL. |
deletion.length |
Integer. Number of deleted nucleotides of pathogenic mutation if it is a DEL. Parameters deletion.sequence and deletion.length can not both be NULL when pathogenic mutation is DEL. Default is NULL. |
insertion.sequence |
Character. Inserted nucleotides of pathogenic mutation if it is an insertion (INS). Default is NULL. |
flanking.window |
Integer. Length of flanking regions (bp) to pathogenic mutation locus for reversion detection. Default is 100. |
minus.strand |
Logical. TRUE if the gene in question is on the reverse strand. Default is FALSE. |
check.wildtype.reads |
Logical. TRUE if assume wildtype reads mapped to pathogenic mutation are restored to wildtype and alternative reversions will be detected from the wildtype reads. Only used for pathogenic mutation targeted gene editing experiment. Default is FALSE. |
A list containing two tables summarizing the reversion detection: 1. The numbers of different types of reads analysed in the input bam file 2. Reversion mutation table including the following columns: rev_id: Unique ID for reversion event rev_freq: Frequency of reversion event rev_type: Type of reversion event, i.e., complement reversion to pathogenic mutation, replacement reversion of pathogenic mutation, or alternative reversion to pathogenic mutation rev_mut_number: Index of each mutation in a reversion event mut_id: Unique ID for reversion mutation chr: Chromosome mut_start_pos: Start position of reversion mutation mut_type: Type of reversion mutation, i.e., SNV, DEL or INS mut_seq: Sequence changes of mutation, i.e., inserted or deleted sequences for indels, or reference and alternative alleles for SNVs mut_length: Length of mut_seq, 0 for SNV mut_hgvs: HGVS Genomic DNA ID of reversion mutation pathog_mut_hgvs: Original pathogenic reversion mutation dist_to_pathog_mut: Distance between original pathogenic mutation and reversion mutation
{ # To detect reversions for BRCA2 mutation "chr13:g.32338763-32338764delAT" bam.file2 <- system.file('extdata', 'toy_alignments_2.bam', package = 'revert') reversions <- getReversions( bam.file = bam.file2, genome.version = "BSgenome.Hsapiens.UCSC.hg38", chromosome = "chr13", pathog.mut.start = 32338763, pathog.mut.type = "DEL", deletion.sequence = "AT", deletion.length = 2, flanking.window = 100, minus.strand = FALSE ) }
{ # To detect reversions for BRCA2 mutation "chr13:g.32338763-32338764delAT" bam.file2 <- system.file('extdata', 'toy_alignments_2.bam', package = 'revert') reversions <- getReversions( bam.file = bam.file2, genome.version = "BSgenome.Hsapiens.UCSC.hg38", chromosome = "chr13", pathog.mut.start = 32338763, pathog.mut.type = "DEL", deletion.sequence = "AT", deletion.length = 2, flanking.window = 100, minus.strand = FALSE ) }