hicdcdiff | R Documentation |
This function calculates differential interactions for a
set of chromosomes across conditions and replicates. You need to install
DESeq2
from Bioconductor to use this function.
hicdcdiff( input_paths, filter_file, output_path, bin_type = "Bins-uniform", binsize = 5000, granularity = 5000, chrs = NULL, Dmin = 0, Dmax = 2e+06, diagnostics = FALSE, DESeq.save = FALSE, fitType = "local" )
input_paths |
a list with names as condition names and values
as paths to |
filter_file |
path to the text file containing columns chr', startI, and startJ denoting the name of the chromosomes and starting coordinates of 2D interaction bins to be compared across conditions, respectively. |
output_path |
the path to the folder and name prefix you want to
place DESeq-processed matrices (in a .txt file), plots
(if |
bin_type |
'Bins-uniform' if uniformly binned by binsize in bp, or 'Bins-RE-sites' if binned by number of restriction enzyme fragment cutsites! |
binsize |
binsize in bp if bin_type='Bins-uniform' (or number of RE fragments if bin_type='Bins-RE-sites'), e.g., default 5000 |
granularity |
Desired distance granularity to base dispersion parameters
on in bp. For uniformly binned analysis
(i.e., |
chrs |
select a subset of chromosomes' e.g., c('chr21','chr22'). Defaults to all chromosomes (except Y and M) in the filter_file. |
Dmin |
minimum distance (included) to check for significant interactions, defaults to 0. Put Dmin=1 to ignore D=0 bins in calculating normalization factors. |
Dmax |
maximum distance (included) to check for significant interactions, defaults to 2e6 or maximum in the data; whichever is minimum. |
diagnostics |
if TRUE, generates diagnostic plots of the normalization factors, geometric means of such factors by distance bin, as well as MA Plots (see DESeq documentation for details about MA plots). Defaults to FALSE. |
DESeq.save |
if TRUE, saves the DESeq objects for each chromosome
as an .rds file in the |
fitType |
follows fitType in |
paths of a list of three entities.
outputpaths
will have differential bins among those in filter_file.
deseq2paths
will have the DESeq2 object stored as an .rds file.
Available if DESeq.save=TRUE
plotpaths
will have diagnostic plots (e.g., MA, dispersion, PCA)
if diagnostics=TRUE
.
outputdir<-paste0(tempdir(check=TRUE),'/') hicdcdiff(input_paths=list(NSD2=c( system.file("extdata", "GSE131651_NSD2_LOW_arima_example.hic", package = "HiCDCPlus"), system.file("extdata", "GSE131651_NSD2_HIGH_arima_example.hic", package = "HiCDCPlus")), TKO=c(system.file("extdata", "GSE131651_TKOCTCF_new_example.hic", package = "HiCDCPlus"), system.file("extdata", "GSE131651_NTKOCTCF_new_example.hic", package = "HiCDCPlus"))), filter_file=system.file("extdata", "GSE131651_analysis_indices.txt.gz", package = "HiCDCPlus"), chrs='chr22', output_path=outputdir, fitType = 'mean', binsize=50000, diagnostics=FALSE)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.