convertToGRanges: Convert InterMineR retrieved genomic information to objects...

Description Usage Arguments Details Value Author(s) See Also Examples

View source: R/convertToGRanges.R

Description

convertToGRanges constitutes a wrapper function for converting genomic locations and their respective annotations to objects of the GRanges class.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
convertToGRanges(
  dataset,
  seqnames,
  start, 
  end,
  names,
  strand,
  columnsAsMetadata = NULL,
  listAsMetadata = NULL,
  seqnamesInterMineChromosome = TRUE
)

Arguments

dataset

a data.frame containing all the genomic information retrieved by InterMineR. Its columns are used to create the GRanges object.

seqnames

the name of a column from dataset or a character vector of length equal to the number of rows of the dataset. Defines the values assigned as seqnames to the GRanges object.

start

the name of a column from dataset or a vector of length equal to the number of rows of the dataset. Defines the genomic coordinates assigned as start to the ranges argument of the GRanges object.

end

the name of a column from dataset or a vector of length equal to the number of rows of the dataset. Defines the genomic coordinates assigned as end to the ranges argument of the GRanges object.

names

the name of a column from dataset or a character vector of length equal to the number of rows of the dataset. Defines the values assigned as names to the ranges argument of the GRanges object.

strand

the name of a column from dataset or a character vector of length equal to the number of rows of the dataset. Defines the values assigned as strand information to the GRanges object.

columnsAsMetadata

a character vector containing the names of the dataset columns that are passed as metadata in the GRanges object.

listAsMetadata

a list of vectors, each of which has length equal to the number of rows of the dataset. The values of the list are passed as metadata to the GRanges object.

seqnamesInterMineChromosome

a logical value indicating whether the values passed as seqnames are InterMineR chromosome.primaryIdentifiers (e.g. "2R", "3R", "X") or not.

Details

The InterMineR package provides a flexible interface to the InterMine web services, which allow for rapid retrieval of various genomic information.

convertToGRanges function is designed to facilitate the conversion of genomic locations and their respective annotations, in the format by which they are retrieved using InterMineR, to an object of the GRanges class.

Value

An object of the GRanges class containing genomic locations and annotations retrieved by InterMineR queries.

Author(s)

InterMine Team

See Also

GRanges, convertToRangedSummarizedExperiment

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
# get FlyMine
im.fly = initInterMine(listMines()["FlyMine"])

# modify template query for Transcription Factor (TF) Binding sites
qTF_Binding = getTemplateQuery(im.fly,"ChromLocation_TFBindingSiteLocationGeneFactor")

qTF_Binding$where[[4]]$value = "1000000"
qTF_Binding$where[[5]]$value = "20000000"

rTF_Binding = runQuery(im.fly, qTF_Binding)

# assign random values for strand of the genomic location retrieved, in InterMine format
rTF_Binding$gene.strand = sample(c("1", "-1", ""), nrow(rTF_Binding), replace = TRUE)

# convert to GRanges object
test = convertToGRanges(
  dataset = rTF_Binding,
  seqnames = rTF_Binding$TFBindingSite.chromosome.primaryIdentifier,
  names = rTF_Binding$TFBindingSite.factor.name,
  start = rTF_Binding$TFBindingSite.chromosomeLocation.start,
  end = rTF_Binding$TFBindingSite.chromosomeLocation.end,
  strand = "gene.strand",
  columnsAsMetadata = c(
    "TFBindingSite.gene.regulatoryRegions.dataSets.dataSource.name",
    "TFBindingSite.factor.primaryIdentifier"),
  listAsMetadata = list(
    c(factor.primaryIdentifier = rTF_Binding$TFBindingSite.factor.primaryIdentifier)
  )
)

# check results
test

intermine/InterMineR documentation built on Jan. 10, 2022, 11:34 p.m.