Package 'm6Aboost'

Title: m6Aboost
Description: This package can help user to run the m6Aboost model on their own miCLIP2 data. The package includes functions to assign the read counts and get the features to run the m6Aboost model. The miCLIP2 data should be stored in a GRanges object. More details can be found in the vignette.
Authors: You Zhou [aut, cre] , Kathi Zarnack [aut]
Maintainer: You Zhou <[email protected]>
License: Artistic-2.0
Version: 1.11.0
Built: 2024-07-21 05:35:59 UTC
Source: https://github.com/bioc/m6Aboost

Help Index


CtoTAssignment for assigning the truncation read counts to the GRanges object with single nucleotide peaks

Description

An function for assigning the CtoT transition read counts from bigWig.

Usage

CtoTAssignment(object, bw_positive, bw_negative, sampleName = "")

Arguments

object

A GRanges object which should contains all the single nucleotide peaks of miCLIP2 experiment.

bw_positive

A path to the bigWig file of C to T transition read counts at the positive strand that output from the preprocess in the m6Aboost pipeline.

bw_negative

A path to the bigWig file of C to T transition read counts at the negative strand that output from the preprocess in the m6Aboost pipeline.

sampleName

The column name that user would like to use for indicating the name of the sample.

Value

A GRanges object with the truncation read counts.

Author(s)

You Zhou

Examples

if (.Platform$OS.type != "windows") {
    testpath <- system.file("extdata", package = "m6Aboost")
    test <- readRDS(file.path(testpath, "test.rds"))
    ctotBw_p <- file.path(testpath, "C2T_positive.bw")
    ctotBw_n <- file.path(testpath, "C2T_negative.bw")
    test <- CtoTAssignment(test, bw_positive=ctotBw_p, bw_negative=ctotBw_n,
        sampleName = "CtoT_WT1")
}

m6Aboost for identify the m6A peaks from the miCILP2 data

Description

An function for calculating the relative signal strength and extracting all the features that required by the m6Aboost model for each peak.

Usage

m6Aboost(object, genome = "", normalization = TRUE)

Arguments

object

A GRanges object which should contains all the single nucleotide peaks of miCLIP2 experiment.

genome

The name of the BSgenome that you are working with. For example "BSgenome.Mmusculus.UCSC.mm10".

normalization

A logical vector which indicates whether you would like normalize the RSS and C to T reads number to the mean value of the training set of the model. This will help to reduce the false positive rate.

Value

A GRanges object with all the information that is required by the m6Aboost model.

Author(s)

You Zhou

Examples

testpath <- system.file("extdata", package = "m6Aboost")
    test_gff3 <- file.path(testpath, "test_annotation.gff3")
    test <- readRDS(file.path(testpath, "test.rds"))
    test<- preparingData(test, test_gff3, colname_reads="WTmean",
        colname_C2T="CtoTmean")

    ## The input of m6Aboost should be the output from preparingData function
    ## Please make sure that the correct BSgenome package have installed
    ## before running motifProfile. For example,
    ## library("BSgenome.Mmusculus.UCSC.mm10")

    test <- m6Aboost(test, "BSgenome.Mmusculus.UCSC.mm10")

preparingData for the miCILP2 data

Description

A function for calculating the relative signal strength and extract the features for running the m6Aboost.

Usage

preparingData(object, annotation, colname_reads = "", colname_C2T = "")

Arguments

object

A GRanges object which should contain all single- nucleotide peaks from the miCLIP2 experiment.

annotation

A path to the annotation file. The format of the annotation file should be a gff3 file and downloaded from https://www.gencodegenes.org/

colname_reads

The name of the metadata column which contains the mean value of the truncation reads number without C to T transition reads.

colname_C2T

The name of the meta data column which contains the mean value of C to T transition read counts.

Value

A GRanges object with all information that is required for running the m6Aboost model.

Author(s)

You Zhou

Examples

testpath <- system.file("extdata", package = "m6Aboost")
    test_gff3 <- file.path(testpath, "test_annotation.gff3")
    test <- readRDS(file.path(testpath, "test.rds"))
    test<- preparingData(test, test_gff3, colname_reads="WTmean",
        colname_C2T="CtoTmean")

truncationAssignment for assigning the truncation read counts to the GRanges object with single nucleotide peaks

Description

An function for assigning the truncation read counts from bigWig to the GRanges peaks.

Usage

truncationAssignment(object, bw_positive, bw_negative, sampleName = "")

Arguments

object

A GRanges object which should contains all the single nucleotide peaks of miCLIP2 experiment.

bw_positive

A path to the bigWig file of truncation read counts at the positive strand that output from the preprocess in the m6Aboost pipeline.

bw_negative

A path to the bigWig file of truncation read counts at the negative strand that output from the preprocess in the m6Aboost pipeline.

sampleName

The column name that user would like to use for indicating the name of the sample.

Value

A GRanges object with the truncation read counts.

Author(s)

You Zhou

Examples

if (.Platform$OS.type != "windows") {
    testpath <- system.file("extdata", package = "m6Aboost")
    test <- readRDS(file.path(testpath, "test.rds"))
    truncationBw_p <- file.path(testpath, "truncation_positive.bw")
    truncationBw_n <- file.path(testpath, "truncation_negative.bw")
    test <- truncationAssignment(test, bw_positive=truncationBw_p,
        bw_negative=truncationBw_n, sampleName = "WT1")

}