Title: | msgbsR: methylation sensitive genotyping by sequencing (MS-GBS) R functions |
---|---|
Description: | Pipeline for the anaysis of a MS-GBS experiment. |
Authors: | Benjamin Mayne |
Maintainer: | Benjamin Mayne <[email protected]> |
License: | GPL-2 |
Version: | 1.31.0 |
Built: | 2024-11-18 04:45:47 UTC |
Source: | https://github.com/bioc/msgbsR |
Determines the sequence around a cut site using a fasta file or BSgenome
checkCuts(cutSites, genome, fasta = FALSE, seq)
checkCuts(cutSites, genome, fasta = FALSE, seq)
cutSites |
A GRanges object containing the locations of the cut sites to be checked for sequence match. The names of the correct cut sites will be returned as a GRanges object. |
genome |
The path to a fasta file or a BSgenome object to check for genomic sequences. |
fasta |
TRUE if a fasta file has been supplied. Default = FALSE |
seq |
The desired recognition sequence that the enzyme should have cut. |
A GRanges object containing the names of the sites that had the correct sequence.
Benjamin Mayne
library(GenomicRanges) library(SummarizedExperiment) library(BSgenome.Rnorvegicus.UCSC.rn6) # Load the positions of possible MspI cut sites data(ratdata) # Extract the cut sites cutSites <- rowRanges(ratdata) # Adjust the cut sites to overlap recognition site on each strand start(cutSites) <- ifelse(test = strand(cutSites) == '+', yes = start(cutSites) - 1, no = start(cutSites) - 2) end(cutSites) <- ifelse(test = strand(cutSites) == '+', yes = end(cutSites) + 2, no = end(cutSites) + 1) correctCuts <- checkCuts(cutSites = cutSites, genome = "rn6", seq = "CCGG")
library(GenomicRanges) library(SummarizedExperiment) library(BSgenome.Rnorvegicus.UCSC.rn6) # Load the positions of possible MspI cut sites data(ratdata) # Extract the cut sites cutSites <- rowRanges(ratdata) # Adjust the cut sites to overlap recognition site on each strand start(cutSites) <- ifelse(test = strand(cutSites) == '+', yes = start(cutSites) - 1, no = start(cutSites) - 2) end(cutSites) <- ifelse(test = strand(cutSites) == '+', yes = end(cutSites) + 2, no = end(cutSites) + 1) correctCuts <- checkCuts(cutSites = cutSites, genome = "rn6", seq = "CCGG")
The GRanges object was created from a list of differentially methylated cut sites from a MS-GBS experiment between two groups of rats that were fed either a control diet or a high fat diet.
data(cuts)
data(cuts)
A GRanges object of length 10.
Positions of MspI cut sites differentially methylated in the prostate on chromosome 20 in Rats.
The data set contains 10 differentially methylated sites in the prostate between rats fed a control or high fat diet.
A GRanges object of length 10.
Determines differential methylated sites from a RangedSummarizedExperiment
diffMeth(se, cateogory, condition1, condition2, block = NULL, cpmThreshold, thresholdSamples)
diffMeth(se, cateogory, condition1, condition2, block = NULL, cpmThreshold, thresholdSamples)
se |
A RangedSummarizedExperiment containing meta data of the samples. |
cateogory |
The heading name in the sample data to be tested for differential methylation. |
condition1 |
The reference group within the cateogory. |
condition2 |
The experimental group within the cateogory. |
block |
The heading name in the sample data if differential methylation is to be tested with a blocking factor. Default is NULL. |
cpmThreshold |
Counts per million threshold of read counts to be filtered out of the analysis. |
thresholdSamples |
Minimum number of samples to contain the counts per million threshold. |
A data frame containing which cut sites that are differenitally methylated.
Benjamin Mayne
# Load data data(ratdata2) top <- diffMeth(se = ratdata2, cateogory = "Group", condition1 = "Control", condition2 = "Experimental", cpmThreshold = 1, thresholdSamples = 1)
# Load data data(ratdata2) top <- diffMeth(se = ratdata2, cateogory = "Group", condition1 = "Control", condition2 = "Experimental", cpmThreshold = 1, thresholdSamples = 1)
Plot a circos representing the cut site locations
plotCircos(cutSites, seqlengths, cutSite.colour, seqlengths.colour)
plotCircos(cutSites, seqlengths, cutSite.colour, seqlengths.colour)
cutSites |
A GRanges object containing the locations of the cut sites to be plotted. |
seqlengths |
An integer with the lengths of the chromosomes. |
cutSite.colour |
The colour of the cut sites. |
seqlengths.colour |
The colour of the chromosomes |
A circos plot showing the locations of the cut sites.
Benjamin Mayne
# load example cut site positions data(cuts) # Obtain the length of chromosome 20 in rn6 library(BSgenome.Rnorvegicus.UCSC.rn6) chr20 <- seqlengths(BSgenome.Rnorvegicus.UCSC.rn6)["chr20"] plotCircos(cutSites = cuts, seqlengths = chr20, cutSite.colour = "red", seqlengths.colour = "blue")
# load example cut site positions data(cuts) # Obtain the length of chromosome 20 in rn6 library(BSgenome.Rnorvegicus.UCSC.rn6) chr20 <- seqlengths(BSgenome.Rnorvegicus.UCSC.rn6)["chr20"] plotCircos(cutSites = cuts, seqlengths = chr20, cutSite.colour = "red", seqlengths.colour = "blue")
Plots the total number of reads vs total number of cut sites per sample
plotCounts(se, cateogory)
plotCounts(se, cateogory)
se |
A RangedSummarizedExperiment containing meta data of the samples. |
cateogory |
The heading name in the sample data to distinguish groups. |
Produces a plot showing the total number reads vs total number of cut sites per sample.
Benjamin Mayne
data(ratdata2) plotCounts(se = ratdata2, cateogory = "Group")
data(ratdata2) plotCounts(se = ratdata2, cateogory = "Group")
A RangedSummarizedExperiment containing read counts generated from a MS-GBS experiment using the restriction enzyme MspI, focusing on chromosome 20 of Rat.
data(ratdata)
data(ratdata)
RangedSummarizedExperiment
ratdata A RangedSummarizedExperiment with 16047 potential MspI cut sites on chromosome 20 in Rat and six samples (3 Control and 3 Experimental).
This dataset contains six prostate samples from rats: 3 control and 3 experimental high fat diet.
RangedSummarizedExperiment
A RangedSummarizedExperiment containing read counts generated from a MS-GBS experiment using the restriction enzyme MspI, focusing on chromosome 20 of Rat. The sites have been checked for the correct recognition site.
data(ratdata2)
data(ratdata2)
RangedSummarizedExperiment
ratdata2 A RangedSummarizedExperiment containing data for 13983 MspI cut sites on chromosome 20 in Rat and six samples (3 Control and 3 Experimental).
This dataset contains six prostate samples from rats: 3 control and 3 experimental high fat diet. The data can be used for differential methylation analyses.
RangedSummarizedExperiment
Imports the raw read counts from sorted and indexed bam file(s)
rawCounts(bamFilepath, threads = 1)
rawCounts(bamFilepath, threads = 1)
bamFilepath |
The path to the location of the bam file(s). |
threads |
The total number of usable threads to be used. Default is 1. |
Produces a RangedSummarizedExperiment. Columns are samples and the rows are cut sites. The cut site IDs are in the format chr:position-position:strand.
Benjamin Mayne, Sam Buckberry
my_path <- system.file("extdata", package = "msgbsR") my_data <- rawCounts(bamFilepath = my_path)
my_path <- system.file("extdata", package = "msgbsR") my_data <- rawCounts(bamFilepath = my_path)