| AddToModBasePropDistFromFeaturePlot | R Documentation |
Add an additional axis with data on a plot generated by DrawModBasePropDistFromFeature function.
AddToModBasePropDistFromFeaturePlot( dPosCountsDistFeatureStart, dPosCountsDistFeatureEnd = NULL, cSubtitleContent, cParamYLabel, cParamColor = "cyan3", cParamType = "l", cParamLwd = 3, cParamLty = 3, lAddAxisOnLeftSide = TRUE, nYLimits = NULL )
dPosCountsDistFeatureStart |
A dataframe containing the data counts by distance to feature positions. |
dPosCountsDistFeatureEnd |
A second dataframe containing the data counts by distance to the second feature positions. If NULL, only 1 position will be plotted (featureStart). Defaults to NULL. |
cSubtitleContent |
A character vector giving the content of the subtitle to add below the title of the plot. |
cParamYLabel |
A character vector giving the content of the label on the new Y axis to add on the plot. |
cParamColor |
The color of the line and new axis to be added. Defaults to "cyan3". |
cParamType |
The type of plot to be added. See "type" argument plot base function. Defaults to "l". |
cParamLwd |
The width of the line/points to be added. Defaults to 3. |
cParamLty |
If a line is plotted, change the type of the line to be added. Defaults to 3. |
lAddAxisOnLeftSide |
If TRUE, add the new axis on the left side of the plot. If FALSE, add the new axis on the right side of the plot. Defaults to TRUE. |
nYLimits |
Numeric vector giving the limits for the new Y axis. Defaults to NULL (will use the minimum and the maximum of both data provided (dPosCountsDistFeatureStart and dPosCountsDistFeatureEnd)). |
# loading genome
myGenome <- Biostrings::readDNAStringSet(system.file(
package = "DNAModAnnot", "extdata",
"ptetraurelia_mac_51_sca171819.fa"
))
# loading annotation
myAnnotations <- rtracklayer::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)
)
# Retrieve, in a list, dataframes of ModBase counts per Distance values from feature positions
myModDistCountsList <- GetDistFromFeaturePos(
grangesAnnotations = myAnnotations,
cSelectFeature = "gene",
grangesData = myGrangesPacBioGFF,
lGetGRangesInsteadOfListCounts = FALSE,
lGetPropInsteadOfCounts = TRUE,
cWhichStrandVsFeaturePos = "both", nWindowSizeAroundFeaturePos = 600,
lAddCorrectedDistFrom5pTo3p = TRUE,
cFeaturePosNames = c("TSS", "TTS")
)
myBaseDistCountsList <- GetDistFromFeaturePos(
grangesAnnotations = myAnnotations,
cSelectFeature = "gene",
grangesData = myGposPacBioCSV,
lGetGRangesInsteadOfListCounts = FALSE,
lGetPropInsteadOfCounts = TRUE,
cWhichStrandVsFeaturePos = "both", nWindowSizeAroundFeaturePos = 600,
lAddCorrectedDistFrom5pTo3p = TRUE,
cFeaturePosNames = c("TSS", "TTS")
)
DrawModBasePropDistFromFeature(
listModCountsDistDataframe = myModDistCountsList,
listBaseCountsDistDataframe = myBaseDistCountsList,
cFeaturePosNames = c("TSS", "TTS"),
cBaseMotif = "A",
cModMotif = "6mA"
)
# Loading Bam data
myBamfile <- Rsamtools::BamFile(file = system.file(
package = "DNAModAnnot", "extdata",
"PTET_MonoNuc_3-2new.pe.sca171819.sorted.bam"
))
myBam_GRanges <- as(GenomicAlignments::readGAlignments(myBamfile), "GRanges")
myBam_GRanges <- GetGposCenterFromGRanges(grangesData = myBam_GRanges)
# Retrieve dataframes of Read center counts per Distance values in a list
myCountsDist_List_bamfile <- GetDistFromFeaturePos(
grangesAnnotations = myAnnotations,
cSelectFeature = "gene",
lGetGRangesInsteadOfListCounts = FALSE,
lGetPropInsteadOfCounts = FALSE,
grangesData = myBam_GRanges,
cWhichStrandVsFeaturePos = "both",
nWindowSizeAroundFeaturePos = 600,
lAddCorrectedDistFrom5pTo3p = TRUE,
cFeaturePosNames = c("TSS", "TTS")
)
# adding new axis to plot from DrawModBasePropDistFromFeature function
AddToModBasePropDistFromFeaturePlot(
dPosCountsDistFeatureStart = myCountsDist_List_bamfile[[1]],
dPosCountsDistFeatureEnd = myCountsDist_List_bamfile[[2]],
cSubtitleContent = "Along with nucleosome center distance (MonoNuc_3-2newreplicate)",
cParamYLabel = "Nucleosome center count (MonoNuc_3-2newreplicate)",
cParamColor = "cyan3",
lAddAxisOnLeftSide = TRUE
)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.