Circular Chromosome Conformation Capture Sequencing (4C-seq) is a sequencing technique enabling the identification of chromatin interactions. Currently, there are several tools available for 4C-seq analysis:
A benchmarking of those tools revealed that none of the tools performed adequately in all evaluated tasks (Walter et al. 2019). To address this we developed fourSynergy, a ensemble algorithm that leverages synergies between the base tools r3cseq, fourSig, peakC, and R.4cker. To this end, all of these tools need to be run, which can be achieved using the fourSynergy pipeline. This pipeline, written in snakemake, includes quality control, interaction calling, and post-processing steps to prepare the output of the base tools for the ensemble algorithm.
If you want to quickly start with testing the package you can jump to the section R analysis and run the fourSynergy R package on testdata.
We recommend the pipeline to be executed in a docker.
To run fourSynergy, you need R (Version 4.5 or higher).
A basic workflow of 4C-seq data is demonstrated in the following on an exemplary dataset. First, we need to run the fourSynergy pipeline.
If you want to quickly start with testing the package you can jump to the section R analysis and run the fourSynergy R package on testdata.
Before starting the R analysis you have to run the fourSynergy snakemake pipeline, which is available at (https://github.com/sophiewind/fourSynergy_pip).
Input for the fourSynergy pipeline are…
4C-seq FASTQ files
config.yaml
Reference genome
The config file is crucial to run fourSynergy. All relevant
experimental information is stored here. You can find an exemplary
info.yaml in inst/extdata/Datasets/Demo/info.yaml.
The pipeline creates a results directory in the directory were its
executed. In the results directory is a directory called like the
author provided in the config file. In this directory, all
results of the base tools and fourSynergy are written.
├── Datasets
├── ...
└── results
└── author
├── alignment -> alignment results
├── basic4Cseq -> Basic4Cseq results
├── foursig -> fourSig results
├── peakC -> peakC results
├── qc -> FASTQC results
├── r3cseq -> r3Cseq results
├── r4cker -> R.4Cker results
└── sia -> standardized tools results
After the pipeline is executed we can start with the analysis in R by loading the fourSynergy R package.
Than the pipeline results can be loaded into an fourSynergy
object. In the provided example, config and the pathes are
stored in extdata. However, in a real-world scenario,
res_path would typically be located in
./results/[DatasetName] and the tracks in
./results/[DatasetName]/alignment.
# Load config
config <- system.file("extdata", "Datasets", "Demo", "info.yaml",
package = "fourSynergy"
)
# Load results path
res_path <- system.file("extdata", "results", "Demo",
package = "fourSynergy"
)
# Load path of .bedGraphs
tracks <- system.file("extdata", "results", "Demo", "alignment",
package = "fourSynergy"
)
sia <- createIa(
res_path = res_path,
config = config,
tracks = tracks
)The fourSynergy object stores information about the single tool calls in condition and in control as well as metadata.
| Slot | Data |
|---|---|
| metadata | Metadata from config |
| expInteractions | Base tool results - experiment |
| ctrlInteractions | Base tool results - control |
| expConsensus | Ensemble results - experiment |
| ctrlConsensus | Ensemble results - control |
| vp | Viewpoint position |
| vfl | Virtual fragment library (VFL) |
| tracks | Path of .bedGraphs |
| differential | Differential ensemble results |
The base tool results are stored in the fourSynergy object in the
slots expInteractions and
ctrlInteractions.
slotNames(sia)
#> [1] "metadata" "expInteractions" "ctrlInteractions" "expConsensus"
#> [5] "ctrlConsensus" "vp" "vfl" "tracks"
#> [9] "differential" "dds"They can be visualized using plotIaIndiviualTools().
There is also the object to pass genes of interest, which are than shown
in the karyoplot.
Further, the base tool results can be visualized in higher resolution
using plotBaseTracks().
#>
#> [[2]]
To perform ensemble calling you need to specify if you want to use the F1 or the AUPRC model.
The ensemble results can than be visualized.
The consensus results can be visualized in higher resolution using
plotConsensusTracks().
#>
#> [[2]]
fourSynergy provides differential interaction analysis to
compare interactions across conditions. The differential interaction
analysis is based on DESeq2 using count data to model the variability in
the data using a negative binomial distribution. In the example dataset
the fitType was set to ‘mean’, by default this value is
‘local’.
(The plots do not look representative here because the BAM file size has been drastically reduced for this example.)
sia <- differentialAnalysis(ia = sia, fitType = "mean")
#> Processing /tmp/RtmpOkCIed/Rinst1aa379057be5/fourSynergy/extdata/results/Demo/alignment/TGFb_1_sorted.bam: I gEt gOoSeBuMPs wHen I seE yOu pilEuPpiNg...
#> Processing /tmp/RtmpOkCIed/Rinst1aa379057be5/fourSynergy/extdata/results/Demo/alignment/TGFb_2_sorted.bam: yOu cAn'T pIlE-Up FaStEr!!!
#> Processing /tmp/RtmpOkCIed/Rinst1aa379057be5/fourSynergy/extdata/results/Demo/alignment/shc_1_sorted.bam: yOu cAn'T pIlE-Up FaStEr!!!
#> Processing /tmp/RtmpOkCIed/Rinst1aa379057be5/fourSynergy/extdata/results/Demo/alignment/shc_2_sorted.bam: yOu cAn'T pIlE-Up FaStEr!!!
#> Warning in DESeqDataSet(se, design = design, ignoreRank): some variables in
#> design formula are characters, converting to factors
#> gene-wise dispersion estimates
#> mean-dispersion relationship
#> final dispersion estimatesThe results of the differential analysis are stored in the
@differential slot of the sia object and can
be visualized using plotDiffIa().
sessionInfo()
#> R version 4.5.2 (2025-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.4 LTS
#>
#> Matrix products: default
#> BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
#> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
#> [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> time zone: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] org.Mm.eg.db_3.22.0 AnnotationDbi_1.73.0 Biobase_2.71.0
#> [4] fourSynergy_0.99.12 GenomicRanges_1.63.1 Seqinfo_1.1.0
#> [7] IRanges_2.45.0 S4Vectors_0.49.0 BiocGenerics_0.57.0
#> [10] generics_0.1.4 BiocStyle_2.39.0
#>
#> loaded via a namespace (and not attached):
#> [1] DBI_1.3.0
#> [2] bitops_1.0-9
#> [3] gridExtra_2.3
#> [4] rlang_1.1.7
#> [5] magrittr_2.0.4
#> [6] biovizBase_1.59.0
#> [7] matrixStats_1.5.0
#> [8] compiler_4.5.2
#> [9] RSQLite_2.4.6
#> [10] GenomicFeatures_1.63.1
#> [11] png_0.1-9
#> [12] vctrs_0.7.1
#> [13] reshape2_1.4.5
#> [14] ProtGenerics_1.43.0
#> [15] stringr_1.6.0
#> [16] pkgconfig_2.0.3
#> [17] crayon_1.5.3
#> [18] fastmap_1.2.0
#> [19] backports_1.5.0
#> [20] XVector_0.51.0
#> [21] labeling_0.4.3
#> [22] Rsamtools_2.27.1
#> [23] rmarkdown_2.30
#> [24] UCSC.utils_1.7.1
#> [25] purrr_1.2.1
#> [26] bit_4.6.0
#> [27] xfun_0.56
#> [28] cachem_1.1.0
#> [29] cigarillo_1.1.0
#> [30] GenomeInfoDb_1.47.2
#> [31] jsonlite_2.0.0
#> [32] blob_1.3.0
#> [33] DelayedArray_0.37.0
#> [34] BiocParallel_1.45.0
#> [35] cluster_2.1.8.2
#> [36] parallel_4.5.2
#> [37] VariantAnnotation_1.57.1
#> [38] R6_2.6.1
#> [39] TxDb.Mmusculus.UCSC.mm10.knownGene_3.10.0
#> [40] bezier_1.1.2
#> [41] bslib_0.10.0
#> [42] stringi_1.8.7
#> [43] RColorBrewer_1.1-3
#> [44] rtracklayer_1.71.3
#> [45] rpart_4.1.24
#> [46] jquerylib_0.1.4
#> [47] Rcpp_1.1.1
#> [48] SummarizedExperiment_1.41.1
#> [49] knitr_1.51
#> [50] base64enc_0.1-6
#> [51] nnet_7.3-20
#> [52] Matrix_1.7-4
#> [53] tidyselect_1.2.1
#> [54] rstudioapi_0.18.0
#> [55] dichromat_2.0-0.1
#> [56] abind_1.4-8
#> [57] yaml_2.3.12
#> [58] codetools_0.2-20
#> [59] curl_7.0.0
#> [60] lattice_0.22-9
#> [61] tibble_3.3.1
#> [62] plyr_1.8.9
#> [63] regioneR_1.43.0
#> [64] withr_3.0.2
#> [65] KEGGREST_1.51.1
#> [66] S7_0.2.1
#> [67] evaluate_1.0.5
#> [68] foreign_0.8-91
#> [69] Biostrings_2.79.5
#> [70] karyoploteR_1.37.0
#> [71] pillar_1.11.1
#> [72] BiocManager_1.30.27
#> [73] MatrixGenerics_1.23.0
#> [74] checkmate_2.3.4
#> [75] TxDb.Hsapiens.UCSC.hg19.knownGene_3.22.1
#> [76] RCurl_1.98-1.17
#> [77] ensembldb_2.35.0
#> [78] ggplot2_4.0.2
#> [79] scales_1.4.0
#> [80] glue_1.8.0
#> [81] lazyeval_0.2.2
#> [82] Hmisc_5.2-5
#> [83] maketools_1.3.2
#> [84] tools_4.5.2
#> [85] BiocIO_1.21.0
#> [86] data.table_1.18.2.1
#> [87] sys_3.4.3
#> [88] BSgenome_1.79.1
#> [89] locfit_1.5-9.12
#> [90] GenomicAlignments_1.47.0
#> [91] buildtools_1.0.0
#> [92] XML_3.99-0.22
#> [93] grid_4.5.2
#> [94] tidyr_1.3.2
#> [95] colorspace_2.1-2
#> [96] htmlTable_2.4.3
#> [97] restfulr_0.0.16
#> [98] Formula_1.2-5
#> [99] cli_3.6.5
#> [100] S4Arrays_1.11.1
#> [101] dplyr_1.2.0
#> [102] AnnotationFilter_1.35.0
#> [103] gtable_0.3.6
#> [104] DESeq2_1.51.7
#> [105] sass_0.4.10
#> [106] digest_0.6.39
#> [107] SparseArray_1.11.11
#> [108] htmlwidgets_1.6.4
#> [109] org.Hs.eg.db_3.22.0
#> [110] rjson_0.2.23
#> [111] farver_2.1.2
#> [112] memoise_2.0.1
#> [113] htmltools_0.5.9
#> [114] lifecycle_1.0.5
#> [115] httr_1.4.8
#> [116] bit64_4.6.0-1
#> [117] bamsignals_1.43.1