The function appends the information about the blocks into the Population Reference SNV Annotation GDS file. The information is extracted from the Population Reference GDS file.

addGeneBlockRefAnnot(
  fileReferenceGDS,
  gdsRefAnnotFile,
  winSize = 10000,
  ensDb,
  suffixBlockName
)

Arguments

fileReferenceGDS

a character string representing the file name of the Reference GDS file. The file must exist.

gdsRefAnnotFile

a character string representing the file name corresponding the Reference SNV Annotation GDS file. The function will open it in write mode and close it after. The file must exist.

winSize

a single positive integer representing the size of the window to use to group the SNVs when the SNVs are in a non-coding region. Default: 10000L.

ensDb

An object with the ensembl genome annotation Default: EnsDb.Hsapiens.v86.

suffixBlockName

a character string that identify the source of the block and that will be added to the block description into the Reference SNV Annotation GDS file, as example: Ensembl.Hsapiens.v86.

Value

The integer OL when the function is successful.

Author

Pascal Belleau, Astrid Deschênes and Alexander Krasnitz

Examples


## Path to the demo pedigree file is located in this package
dataDir <- system.file("extdata", package="RAIDS")

fileAnnotGDS <- file.path(tempdir(), "ex1_good_small_1KG_Ann_GDS.gds")

## Required library
if (requireNamespace("EnsDb.Hsapiens.v86", quietly=TRUE)) {

    file.copy(file.path(dataDir, "tests",
        "ex1_NoBlockGene.1KG_Annot_GDS.gds"), fileAnnotGDS)

    ## Making a "short cut" on the ensDb object
    edb <- EnsDb.Hsapiens.v86::EnsDb.Hsapiens.v86

    ## GDS Reference file
    fileReferenceGDS  <- file.path(dataDir, "tests",
                "ex1_good_small_1KG.gds")

    # \donttest{


        ## Append information associated to blocks
        addGeneBlockRefAnnot(fileReferenceGDS=fileReferenceGDS,
            gdsRefAnnotFile=fileAnnotGDS,
            ensDb=edb,
            suffixBlockName="EnsDb.Hsapiens.v86")

        gdsAnnot1KG <- openfn.gds(fileAnnotGDS)
        print(gdsAnnot1KG)
        print(read.gdsn(index.gdsn(gdsAnnot1KG, "block.annot")))

       closefn.gds(gdsAnnot1KG)
    # }

    ## Remove temporary file
    unlink(fileAnnotGDS, force=TRUE)

}
#> File: /tmp/RtmplgNUHa/ex1_good_small_1KG_Ann_GDS.gds (392.5K)
#> +    [  ]
#> |--+ phase   { Bit2 11000x156 LZ4_ra(21.5%), 90.2K }
#> |--+ block.annot   [ data.frame ] *
#> |  |--+ block.id   { Str8 7, 119B }
#> |  \--+ block.desc   { Str8 7, 384B }
#> \--+ block   { Int32 11000x7 LZ4_ra(22.5%), 67.6K }
#>                   block.id
#> 1            EAS.0.05.500k
#> 2            EUR.0.05.500k
#> 3            AFR.0.05.500k
#> 4            AMR.0.05.500k
#> 5            SAS.0.05.500k
#> 6  Gene.EnsDb.Hsapiens.v86
#> 7 GeneS.EnsDb.Hsapiens.v86
#>                                                      block.desc
#> 1         EAS populationblock base on SNP 0.05 and windows 500k
#> 2         EUR populationblock base on SNP 0.05 and windows 500k
#> 3         AFR populationblock base on SNP 0.05 and windows 500k
#> 4         AMR populationblock base on SNP 0.05 and windows 500k
#> 5         SAS populationblock base on SNP 0.05 and windows 500k
#> 6 List of blocks including overlapping genes EnsDb.Hsapiens.v86
#> 7           List of blocks of split by genes EnsDb.Hsapiens.v86