Skip to contents

Export the predicted gene-pair interactions calculated by upstream function linkGenesAndPeaks into an Interact Track file which is compatible with UCSC Genome Browser.

Usage

exportInteractTrack(
  corrMat,
  pathToCoords,
  useGenes = NULL,
  outputPath = getwd()
)

Arguments

corrMat

A sparse matrix of correlation with peak names as rows and gene names as columns.

pathToCoords

Path to the gene coordinates file.

useGenes

Character vector of gene names to be exported. Default NULL uses all genes available in corrMat.

outputPath

Path of filename where the output file will be stored. If a folder, a file named "Interact_Track.bed" will be created. Default current working directory.

Value

No return value. A file located at outputPath will be created.

Examples

# \donttest{
bmmc <- normalize(bmmc)
#>  Normalizing datasets "rna"
#>  Normalizing datasets "atac"
#>  Normalizing datasets "atac" ... done
#> 
#>  Normalizing datasets "rna"

#>  Normalizing datasets "rna" ... done
#> 
bmmc <- selectGenes(bmmc)
#>  Selecting variable features for dataset "rna"
#>  ... 83 features selected out of 172 shared features.
#>  Selecting variable features for dataset "atac"
#>  ... 126 features selected out of 172 shared features.
#>  Finally 135 shared variable features are selected.
bmmc <- scaleNotCenter(bmmc)
#>  Scaling dataset "rna"
#>  Scaling dataset "rna" ... done
#> 
#>  Scaling dataset "atac"
#>  Scaling dataset "atac" ... done
#> 
if (requireNamespace("RcppPlanc", quietly = TRUE) &&
    requireNamespace("GenomicRanges", quietly = TRUE) &&
    requireNamespace("IRanges", quietly = TRUE) &&
    requireNamespace("psych", quietly = TRUE)) {
    bmmc <- runINMF(bmmc)
    bmmc <- quantileNorm(bmmc)
    bmmc <- normalizePeak(bmmc)
    bmmc <- imputeKNN(bmmc, reference = "atac", queries = "rna")
    corr <- linkGenesAndPeaks(
        bmmc, useDataset = "rna",
        pathToCoords = system.file("extdata/hg19_genes.bed", package = "rliger")
    )
    resultPath <- tempfile()
    exportInteractTrack(
        corrMat = corr,
        pathToCoords = system.file("extdata/hg19_genes.bed", package = "rliger"),
        outputPath = resultPath
    )
    head(read.table(resultPath, skip = 1))
}
#>  Using largest dataset of recommended type as reference: "rna" with 340 cells
#>  Normalizing peak of dataset: "atac"
#>  Normalizing peak of dataset: "atac" ... done
#> 
#>  Imputing 1 query dataset: "rna"
#>  from reference dataset: "atac"
#>  Normalizing peak of dataset: "rna"
#>  Normalizing peak of dataset: "rna" ... done
#> 
#>  172 genes to be tested against 995 peaks
#>  Calculating correlation for gene-peak pairs...
#>  ■■■■■■■■■■■■■■■■■■■■■■■■■■■■■     93% | ETA:  0s
#>  ■■■■■■■■■■■■■■■■■■■■■■■■■■■■■■■  100% | ETA:  0s
#> ! Totally 141 selected genes do not have significant correlated peaks, out of
#> 172 selected genes
#>  Result written at: /private/var/folders/k9/nwtr_c_57kd43cmbtz82ywlr0000gn/T/RtmpAnYSmu/file31412559055c
#>     V1        V2        V3                             V4 V5         V6 V7 V8
#> 1 chr1 155928566 155929066  LMNA/chr1:155928566-155929066  0  0.2997087  .  5
#> 2 chr1 220868226 220868726 MARC2/chr1:220868226-220868726  0  0.1840661  .  5
#> 3 chr2  20382204  20382704    SDC1/chr2:20382204-20382704  0  0.5781944  .  5
#> 4 chr2  20412152  20412652    SDC1/chr2:20412152-20412652  0  0.1810287  .  5
#> 5 chr2 232223879 232224379 ITM2C/chr2:232223879-232224379  0 -0.1559540  .  5
#> 6 chr3  52090371  52090871   STAB1/chr3:52090371-52090871  0  0.1717412  .  5
#>     V9       V10       V11 V12 V13  V14       V15       V16   V17 V18
#> 1 chr1 155928566 155928567   .   . chr1 156052368 156052369  LMNA   .
#> 2 chr1 220868226 220868227   .   . chr1 220921675 220921676 MARC2   .
#> 3 chr2  20382204  20382205   .   . chr2  20400557  20400558  SDC1   .
#> 4 chr2  20412152  20412153   .   . chr2  20400557  20400558  SDC1   .
#> 5 chr2 232223879 232223880   .   . chr2 231729620 231729621 ITM2C   .
#> 6 chr3  52090371  52090372   .   . chr3  52529355  52529356 STAB1   .
# }