Title: | Fast Implementation of Reading/Dump for .hic Files |
---|---|
Description: | API for fast data extraction for .hic files that provides programmatic access to the matrices. It doesn't store the pointer data for all the matrices, only the one queried, and currently we are only supporting matrices (not vectors). |
Authors: | Neva Cherniavsky Durand [aut, cre], Muhammad Saad Shamim [aut], Aiden Lab [cph] |
Maintainer: | Neva Cherniavsky Durand <[email protected]> |
License: | MIT + file LICENSE |
Version: | 0.0.9 |
Built: | 2024-11-12 05:26:24 UTC |
Source: | https://github.com/aidenlab/straw |
Function for reading basepair resolutions from .hic file
readHicBpResolutions(fname)
readHicBpResolutions(fname)
fname |
path to .hic file |
Vector of basepair resolutions
readHicBpResolutions(system.file("extdata", "test.hic", package = "strawr"))
readHicBpResolutions(system.file("extdata", "test.hic", package = "strawr"))
Function for reading chromosomes from .hic file
readHicChroms(fname)
readHicChroms(fname)
fname |
path to .hic file |
Data frame of chromosome names and lengths
readHicChroms(system.file("extdata", "test.hic", package = "strawr"))
readHicChroms(system.file("extdata", "test.hic", package = "strawr"))
Function for reading available normalizations from .hic file
readHicNormTypes(fname)
readHicNormTypes(fname)
fname |
path to .hic file |
Vector of available normalizations
readHicNormTypes(system.file("extdata", "test.hic", package = "strawr"))
readHicNormTypes(system.file("extdata", "test.hic", package = "strawr"))
fast C++ implementation of dump. Not as fully featured as the Java version. Reads the .hic file, finds the appropriate matrix and slice of data, and outputs as data.frame in sparse upper triangular format. Currently only supporting matrices.
straw(norm, fname, chr1loc, chr2loc, unit, binsize, matrix = "observed")
straw(norm, fname, chr1loc, chr2loc, unit, binsize, matrix = "observed")
norm |
Normalization to apply. Must be one of NONE/VC/VC_SQRT/KR. VC is vanilla coverage, VC_SQRT is square root of vanilla coverage, and KR is Knight-Ruiz or Balanced normalization. |
fname |
path to .hic file |
chr1loc |
first chromosome location |
chr2loc |
second chromosome location |
unit |
BP (BasePair) or FRAG (FRAGment) |
binsize |
The bin size. By default, for BP, this is one of <2500000, 1000000, 500000, 250000, 100000, 50000, 25000, 10000, 5000> and for FRAG this is one of <500, 200, 100, 50, 20, 5, 2, 1>. |
matrix |
Type of matrix to output. Must be one of observed/oe/expected. observed is observed counts, oe is observed/expected counts, expected is expected counts. |
Usage: straw <NONE/VC/VC_SQRT/KR> <hicFile(s)> <chr1>[:x1:x2] <chr2>[:y1:y2] <BP/FRAG> <binsize> [observed/oe/expected]
Data.frame of a sparse matrix of data from hic file. x,y,counts
straw("NONE", system.file("extdata", "test.hic", package = "strawr"), "1", "1", "BP", 2500000)
straw("NONE", system.file("extdata", "test.hic", package = "strawr"), "1", "1", "BP", 2500000)
API for fast data extraction for .hic files that provides programmatic access to the matrices. It doesn't store the pointer data for all the matrices, only the one queried, and currently we are only supporting matrices (not vectors).