searchVariants: searchVariants function

Description Usage Arguments Details Value References See Also Examples

View source: R/searchVariants.R

Description

Search for variants by genomic ranges (lines of VCF files).

Usage

1
2
3
searchVariants(host, variantSetId, referenceName, start, end,
  callSetIds = character(), nrows = Inf, responseSize = NA_integer_,
  asVCF = TRUE)

Arguments

host

URL of GA4GH API data server.

variantSetId

The variant set to search.

referenceName

Required. Only return variants on this reference.

start

Required. The beginning of the window (1-based, inclusive) for which overlapping variants should be returned. Genomic positions are non-negative integers less than reference length. Requests spanning the join of circular genomes are represented as two requests one on each side of the join (position 1).

end

Required. The end of the window (1-based, inclusive) for which overlapping variants should be returned.

callSetIds

Only return variant calls which belong to callsets with these IDs. If unspecified, return all variants and no variant call objects.

nrows

Number of rows of the data frame returned by this function. If not defined, the function will return all entries. If the number of available entries is less than the value of this this parameter, the function will silently return only the available entries.

responseSize

Specifies the number of entries to be returned by the server until reach the number of rows defined in nrows parameter or until get all available entries. If not defined, the server will return the allowed maximum reponse size. Increasing this the value of this parameter will reduce the number of requests and reducing the time required. The will not respect this parameter if the value if larger than its maximum response size.

asVCF

If TRUE the function will return an VCF with header (default), otherwise it will return an DataFrame.

Details

This function maps to POST host/variants/search.

Value

VCF object (when asVCF = TRUE) or DataFrame object (otherwise).

References

Official documentation.

See Also

DataFrame, getVariant, searchVariantsByGRanges, VCF, makeVCFFromGA4GHResponse

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
host <- "http://1kgenomes.ga4gh.org/"
## Not run: 
datasetId <- searchDatasets(host, nrows = 1)$id
variantSetId <- searchVariantSets(host, datasetId, nrows = 1)$id
searchVariants(host, variantSetId, referenceName = "1",
    start = 15000, end = 16000)

searchVariants(host, variantSetId, referenceName = "1",
    start = 15000, end = 16000, asVCF = FALSE)

## End(Not run)

labbcb/GA4GHclient documentation built on May 20, 2019, 7:32 p.m.