GetModBaseCountsWithinFeature | R Documentation |
Cut the Annotation provided into x windows of relative size and return, for each window, the counts (and counts per KiloBase pairs (kbp)) of the base modified (Mod) and the base letter of the modified base (Base). Example: for Mod="6mA", Base="A"; for Mod="5mC", Base="C".
GetModBaseCountsWithinFeature( grangesAnnotations, nWindowsNb = 20, grangesModPos, gposModTargetBasePos, lIgnoreStrand = FALSE )
grangesAnnotations |
A GRanges object containing the annotation for the genome assembly corresponding to the grangesModPos and gposModTargetBasePos provided. The Genomic features categories must be in a column named "type". |
nWindowsNb |
The number of output windows by feature. The annotation provided (with input ranges) will be cut into this number of output windows. Each output window will have the same size as the other output windows from the same input range. Defaults to 20. |
grangesModPos |
A GRanges object containing Modifications Positions data to be counted. |
gposModTargetBasePos |
A GRanges or GPos object containing Base Positions (which can be targeted by the modification) to be counted. |
lIgnoreStrand |
If TRUE, Mod and Base will be counted independently of the strand of each feature. If FALSE, only Mod and Base on the same strand as the feature will be counted. Defaults to FALSE. |
A GRanges object based on grangesAnnotations with the counts:
Modcount: The number of "Mod" within this window of this feature.
ModcountSum: Total number of "Mod" within all windows of this feature.
Modprop: (Modcount / ModcountSum)*100
Basecount: The number of "Base" within this window of this feature.
BasecountSum: Total number of "Base" within all windows of this feature.
Baseprop: (Basecount/BasecountSum)*100
# loading genome myGenome <- Biostrings::readDNAStringSet(system.file( package = "DNAModAnnot", "extdata", "ptetraurelia_mac_51_sca171819.fa" )) # loading annotation library(rtracklayer) myAnnotations <- readGFFAsGRanges(system.file( package = "DNAModAnnot", "extdata", "ptetraurelia_mac_51_annotation_v2.0_sca171819.gff3" )) # Preparing a grangesPacBioGFF and a grangesPacBioCSV datasets myGrangesPacBioGFF <- ImportPacBioGFF( cPacBioGFFPath = system.file( package = "DNAModAnnot", "extdata", "ptetraurelia.modifications.sca171819.gff" ), cNameModToExtract = "m6A", cModNameInOutput = "6mA", cContigToBeAnalyzed = names(myGenome) ) myGposPacBioCSV <- ImportPacBioCSV( cPacBioCSVPath = system.file( package = "DNAModAnnot", "extdata", "ptetraurelia.bases.sca171819.csv" ), cSelectColumnsToExtract = c( "refName", "tpl", "strand", "base", "score", "ipdRatio", "coverage" ), lKeepExtraColumnsInGPos = TRUE, lSortGPos = TRUE, cContigToBeAnalyzed = names(myGenome) ) myGposPacBioCSV <- myGposPacBioCSV[myGposPacBioCSV$base == "A"] # Retrieve annotations with "Mod" and "Base" counts (and counts per kbp) myAnn_ModBase_counts_by_window <- GetModBaseCountsWithinFeature( grangesAnnotations = myAnnotations[myAnnotations$type == "gene", ], grangesModPos = myGrangesPacBioGFF, gposModTargetBasePos = myGposPacBioCSV, nWindowsNb = 20 ) myAnn_ModBase_counts_by_window
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.