This function adapts the chicane
function from the ChICANE
package to work with the linkSet
object format. It runs the full method for
detecting significant interactions in capture Hi-C experiments.
Usage
run_chicane(linkSet, ...)
# S4 method for class 'linkSet'
run_chicane(
linkSet,
replicate.merging.method = "sum",
distribution = "negative-binomial",
include.zeros = "none",
bait.filters = c(0, 1),
target.filters = c(0, 1),
distance.bins = NULL,
multiple.testing.correction = c("bait-level", "global"),
adjustment.terms = NULL,
remove.adjacent = FALSE,
temp.directory = NULL,
keep.files = FALSE,
maxit = 100,
epsilon = 1e-08,
cores = 1,
trace = FALSE,
verbose = FALSE
)
Arguments
- linkSet
A linkSet object containing interaction data
- ...
Additional arguments passed to methods
- replicate.merging.method
Method for merging replicates (default: 'sum')
- distribution
Distribution to use for modeling (default: 'negative-binomial')
- include.zeros
How to handle zero counts (default: 'none')
- bait.filters
Vector of length 2 for bait filtering thresholds (default: c(0,1))
- target.filters
Vector of length 2 for target filtering thresholds (default: c(0,1))
- distance.bins
Number of distance bins (default: NULL)
- multiple.testing.correction
Method for multiple testing correction (default: 'bait-level')
- adjustment.terms
Additional terms for model adjustment (default: NULL)
- remove.adjacent
Whether to remove adjacent fragments (default: FALSE)
- temp.directory
Directory for temporary files (default: NULL)
- keep.files
Whether to keep temporary files (default: FALSE)
- maxit
Maximum iterations for model fitting (default: 100)
- epsilon
Convergence threshold (default: 1e-8)
- cores
Number of CPU cores to use (default: 1)
- trace
Whether to print trace information (default: FALSE)
- verbose
Whether to print progress information (default: FALSE)
Value
A linkSet object with additional columns:
expectedThe expected number of reads linking fragments under the fitted model
p.valueP-value for test of observed vs expected read counts
q.valueFDR-corrected p-value
Examples
# Create example data
gr1 <- GRanges(seqnames = c("chr1", "chr3", "chr3"),
ranges = IRanges(start = c(1000, 2000, 3000), width = 100),
strand = "+", symbol = c("BRCA1", "TP53", "NONEXISTENT"))
gr2 <- GRanges(seqnames = c("chr1", "chr2", "chr3"),
ranges = IRanges(start = c(5000, 6000, 7000), width = 100),
strand = "+")
ls <- linkSet(gr1, gr2, specificCol = "symbol")
# Annotate and prepare data
annotated_ls <- suppressWarnings(
annotatePromoter(ls, genome = "hg38", upstream = 500, overwrite = TRUE)
)
annotated_ls <- countInteractibility(annotated_ls)
annotated_ls <- linkSet::pairdist(annotated_ls)
# Run analysis
result_ls <- run_chicane(
annotated_ls,
replicate.merging.method = 'sum',
bait.filters = c(0, 1),
target.filters = c(0, 1),
distance.bins = NULL,
multiple.testing.correction = 'bait-level',
verbose = TRUE
)
#> [1] "Not found column 'bait.id', adding bait name as default."
#> [1] "Not found column 'bait.to.bait' Set 'FALSE' as default."
#> FITTING MODEL
#> [1] "fitting model...."
#>
#> Attaching package: ‘data.table’
#> The following object is masked from ‘package:SummarizedExperiment’:
#>
#> shift
#> The following object is masked from ‘package:GenomicRanges’:
#>
#> shift
#> The following object is masked from ‘package:IRanges’:
#>
#> shift
#> The following objects are masked from ‘package:S4Vectors’:
#>
#> first, second
#> *
#> Warning: no observations informative at iteration 1
#> Warning: no observations informative at iteration 1
#> Warning: glm.fit: algorithm did not converge
#>
#> trans interactions