DrawModBasePropByFeature: DrawModBasePropByFeature Function (ModAnnot)

View source: R/ModAnnot.R

DrawModBasePropByFeatureR Documentation

DrawModBasePropByFeature Function (ModAnnot)

Description

Return a plot describing, between the features to be displayed, the proportion (or proportion 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".

Usage

DrawModBasePropByFeature(
  grangesAnnotationsWithCounts,
  cFeaturesToCompare = c("gene", "intergenic"),
  lUseCountsPerkbp = FALSE,
  cBaseMotif,
  cModMotif
)

Arguments

grangesAnnotationsWithCounts

A GRanges object based on grangesAnnotations with the counts:

  • Modcount: The number of "Mod" within this feature.

  • Modcount_perkbp: (The number of "Mod" within this feature divided by the size of the feature)*1000.

  • Basecount: The number of "Base" within this feature.

  • Basecount_perkbp: (The number of "Base" within this feature divided by the size of the feature)*1000.

The Genomic features categories must be in a column named "type".

cFeaturesToCompare

Names of the features to be displayed and compared. Defaults to c("gene", "intergenic").

lUseCountsPerkbp

If TRUE, use counts per kbp (Modcount_perkbp and Basecount_perkbp) to calculate the proportion between features displayed. If FALSE, use counts (Modcount and Basecount) to calculate the proportion between features displayed. Defaults to FALSE.

cBaseMotif

The name of the motif with the base letter of the modified base.

cModMotif

The name of the motif with the modification in the output.

Examples

# loading genome
myGenome <- Biostrings::readDNAStringSet(system.file(
  package = "DNAModAnnot", "extdata",
  "ptetraurelia_mac_51_sca171819.fa"
))
myGrangesGenome <- GetGenomeGRanges(myGenome)

# loading annotation
library(rtracklayer)
myAnnotations <- readGFFAsGRanges(system.file(
  package = "DNAModAnnot", "extdata",
  "ptetraurelia_mac_51_annotation_v2.0_sca171819.gff3"
))
myAnnotations <- PredictMissingAnnotation(
  grangesAnnotations = myAnnotations,
  grangesGenome = myGrangesGenome,
  cFeaturesColName = "type",
  cGeneCategories = c("gene"),
  lAddIntronRangesUsingExon = TRUE
)

# 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 <- GetModBaseCountsByFeature(
  grangesAnnotations = myAnnotations,
  grangesModPos = myGrangesPacBioGFF,
  gposModTargetBasePos = myGposPacBioCSV
)

DrawModBasePropByFeature(
  grangesAnnotationsWithCounts = myAnn_ModBase_counts,
  cFeaturesToCompare = c("gene", "intergenic"),
  lUseCountsPerkbp = TRUE,
  cBaseMotif = "A",
  cModMotif = "6mA"
)

AlexisHardy/DNAModAnnot documentation built on Feb. 27, 2023, 12:03 a.m.