SNP Effect Matrices
SEMplR uses SNP Effect Matrices (SEMs) to score potential motifs.
These are matrices contain binding affinity scores and have rows equal
to the length of the motif and a column for each nucleotide. SEMs are
produced by SEMpl, but a default set of 223 are included with this
package in the sc
data object. A full list of the
transcription factors included in this default set can be found here
or by running semData(sc)
.
The default SEM collection can be loaded with:
data(sc)
Printing the SNPEffectMatrixCollection
, we can see that
we have 223 SEMs and 13 meta data features for each.
sc
#> An object of class SNPEffectMatrixCollection
#> sems(223): AP2B_HUMAN.SK-N-SH, ARNT_HUMAN.GM12878 ... ZNF18_HUMAN.HEK293, ZSCAN4_secondary
#> semData(13): transcription_factor, ensembl_id, ... , PWM_source, SEM_KEY
We can access all SEMs in collection the with the function
sems()
or some subset of matrices by specifying a vector of
semIds in the semId
parameter.
sems(sc, semId = c("AP2B_HUMAN.SK-N-SH", "ZSCAN4_secondary"))
#> $`AP2B_HUMAN.SK-N-SH`
#> An object of class SNPEffectMatrix
#> semId: AP2B_HUMAN.SK-N-SH
#> baseline: -1.161048
#> sem:
#> A C G T
#> <num> <num> <num> <num>
#> 1: -1.093659 -0.228215 -0.035353 -0.949757
#> 2: -1.198571 0.019236 -1.057251 -1.145949
#> 3: -1.021307 0.008023 -1.273623 -0.746918
#> 4: -0.808380 -0.281746 -0.849156 0.093359
#> 5: -0.319849 -0.071797 -0.058358 -0.085237
#> 6: 0.119797 -0.640756 -0.336054 -0.787809
#> 7: -0.496324 -1.111405 0.034943 -0.894114
#> 8: -1.069881 -1.178852 0.005380 -1.079589
#> 9: -0.957910 -0.095025 -0.061232 -1.025970
#> 10: -0.147762 -0.046611 -0.308033 -0.113474
#>
#> $ZSCAN4_secondary
#> An object of class SNPEffectMatrix
#> semId: ZSCAN4_secondary
#> baseline: -2.911018
#> sem:
#> A C G T
#> <num> <num> <num> <num>
#> 1: 0.339681 -1.6258400 -0.852570 -1.40307
#> 2: -1.858890 0.2245480 -1.803160 -1.10145
#> 3: 0.259177 -1.8665800 -0.744934 -1.62341
#> 4: -2.077230 0.1970890 -1.971940 -1.21187
#> 5: 0.198878 -2.2059500 -0.760569 -1.76139
#> 6: -2.321170 0.1354730 -2.300640 -1.39237
#> 7: 0.157553 -2.4685400 -0.741112 -1.84194
#> 8: -2.482910 0.0893564 -2.406690 -1.56081
#> 9: 0.136749 -2.4335100 -0.704496 -1.94128
#> 10: -2.352060 0.1257340 -2.211630 -1.60719
#> 11: 0.187489 -2.1923800 -0.641008 -1.76197
#> 12: -2.109030 0.1677660 -1.974360 -1.51917
#> 13: 0.240220 -1.9354700 -0.655348 -1.54698
#> 14: -1.748380 0.2142480 -1.693530 -1.38530
#> 15: 0.289375 -1.6708800 -0.698449 -1.42560
#> 16: -1.483480 0.2650430 -1.541780 -1.31953
We can view the SEM meta data slot with the semData
function.
semData(sc)
#> Key: <SEM_KEY>
#> transcription_factor ensembl_id ebi_complex_ac uniprot_ac
#> <char> <char> <char> <char>
#> 1: TFAP2B ENSG00000008196 Q92481
#> 2: ARNT ENSG00000143437 P27540
#> 3: ATF1 ENSG00000123268 P18846
#> 4: ATF2 ENSG00000115966 P15336
#> 5: ATF3 ENSG00000162772 P18847
#> ---
#> 219: ZBTB7A ENSG00000178951 O95365
#> 220: ZFX ENSG00000005889 P17010
#> 221: ZNF281 ENSG00000162702 Q9Y2X9
#> 222: ZNF18 ENSG00000154957 P17022
#> 223: ZSCAN4 ENSG00000180532 Q8NAM6
#> PWM_id SEM SEM_baseline cell_type
#> <char> <char> <num> <char>
#> 1: AP2B_HUMAN.H11MO.0.B AP2B_HUMAN.SK-N-SH.sem -1.1610478 SK-N-SH
#> 2: ARNT_HUMAN.H11MO.0.B ARNT_HUMAN.GM12878.sem -1.8632604 GM12878
#> 3: ATF1_HUMAN.H11MO.0.B ATF1_HUMAN.K562.sem -2.9718855 K562
#> 4: ATF2_HUMAN.H11MO.0.B ATF2_HUMAN.HepG2.sem -1.4186930 HepG2
#> 5: ATF3_HUMAN.H11MO.0.A ATF3_HUMAN.HepG2.sem -2.1588687 HepG2
#> ---
#> 219: ZBT7A_HUMAN.H11MO.0.A ZBT7A_HUMAN.HepG2.sem -0.4056000 HepG2
#> 220: ZFX_HUMAN.H11MO.0.A ZFX_HUMAN.HepG2.sem -0.2478717 HepG2
#> 221: ZN281_HUMAN.H11MO.0.A ZN281_HUMAN.HepG2.sem -0.3926150 HepG2
#> 222: ZNF18_HUMAN.H11MO.0.C ZNF18_HUMAN.HEK293.sem -1.6459529 HEK293
#> 223: Zscan4_pwm_secondary ZSCAN4_secondary.sem -2.9110180 HEK293
#> neg_log10_pval chip_ENCODE_accession dnase_ENCODE_accession PWM_source
#> <num> <char> <char> <char>
#> 1: 6.476702 ENCFF303ORH ENCFF752OZB HOCOMOCOv11
#> 2: 12.081160 ENCFF538DNI ENCFF759OLD HOCOMOCOv11
#> 3: 56.056680 ENCFF885JPB ENCFF274YGF HOCOMOCOv11
#> 4: 27.821210 ENCFF899JAX ENCFF897NME HOCOMOCOv11
#> 5: 14.618530 ENCFF197JBM ENCFF897NME HOCOMOCOv11
#> ---
#> 219: 3.684286 ENCFF439DRT ENCFF897NME HOCOMOCOv11
#> 220: 3.443155 ENCFF329DCV ENCFF897NME HOCOMOCOv11
#> 221: 3.428189 ENCFF948PYK ENCFF897NME HOCOMOCOv11
#> 222: 3.793244 ENCFF888BHA ENCFF285OXK HOCOMOCOv11
#> 223: 11.973690 ENCFF485WHZ ENCFF910QHN UNIPROBE
#> SEM_KEY
#> <char>
#> 1: AP2B_HUMAN.SK-N-SH
#> 2: ARNT_HUMAN.GM12878
#> 3: ATF1_HUMAN.K562
#> 4: ATF2_HUMAN.HepG2
#> 5: ATF3_HUMAN.HepG2
#> ---
#> 219: ZBT7A_HUMAN.HepG2
#> 220: ZFX_HUMAN.HepG2
#> 221: ZN281_HUMAN.HepG2
#> 222: ZNF18_HUMAN.HEK293
#> 223: ZSCAN4_secondary
Scoring Binding
plotSEM(sc, motif = "AP2B_HUMAN.SK-N-SH",
motifSeq = "GCTTTGAGGC", highlight = 1)
plotSEM(sc, motif = "AP2B_HUMAN.SK-N-SH",
motifSeq = "GCTTTCAGGC", highlight = 1, hcol = "red")
Scoring Variants
Prepare inputs
First, we will define the variants we want to score. Variants can be
supplied to SEMplR as either a VRanges
or a
GRanges
object. If using VRanges
alleles must
be stored in the ref
and alt
parameters. If
using GRanges
the alleles should be stored in seperate
metadata columns.
Insertions and deletions should be specified by an empty character vector in the appropriate allele.
Optionally, an id
column can be defined in the object
metadata to be used as a unique identifier for each variant.
Scoring
The scoreVariants
function scores each allele of each
variant against every SEM provided. From the resulting data object, we
can see that we scored 2 variants, the SEM meta data stored in the
object has 13 fields, and we have 446 rows in the scores
results table.
There is a row in the scores
table for each variant/SEM
combination. Here we scored 2 variants x 223 SEMs to get 446 rows.
sempl_results <- scoreVariants(vr = vr,
semList = sc,
bs_genome_obj=BSgenome.Hsapiens.UCSC.hg38::Hsapiens)
sempl_results
#> An object of class SemplScores
#> variants(2): variant1, variant2
#> semData(13): transcription_factor, ensembl_id, ... , PWM_source, SEM_KEY
#> scores(446):
#> varId semId refSeq altSeq refScore
#> <char> <char> <char> <char> <num>
#> 1: variant1 AP2B_HUMAN.SK-N-SH GTTTAACAAT TCTTTAACAA -5.287583
#> 2: variant1 ARNT_HUMAN.GM12878 TTGTTTAAC TGTTGTTCT -5.048425
#> 3: variant1 ATF1_HUMAN.K562 ACTTGTTGTTG CTTTAACAATA -7.736531
#> 4: variant1 ATF2_HUMAN.HepG2 CTTGTTGTTGT CTTGTTGTTCT -5.049124
#> 5: variant1 ATF3_HUMAN.HepG2 TGTTGTTTAAC TGTTCTTTAAC -9.214165
#> ---
#> 442: variant2 ZBT7A_HUMAN.HepG2 TTTAATCCT TATAATCCT -2.821984
#> 443: variant2 ZFX_HUMAN.HepG2 TTAATCCTCT ATAATCCTCT -1.407772
#> 444: variant2 ZN281_HUMAN.HepG2 CTTGGGCAAATTATT CTTGGGCAAATTATA -5.109413
#> 445: variant2 ZNF18_HUMAN.HEK293 TATTTAATCCTC TATATAATCCTC -6.516791
#> 446: variant2 ZSCAN4_secondary TCTTGGGCAAATTATT TTATATAATCCTCTAA -17.314249
#> altScore refNorm altNorm refVarIndex altVarIndex
#> <num> <num> <num> <int> <int>
#> 1: -4.330368 -0.9427482 -0.8888423 20 19
#> 2: -5.948900 -0.8900564 -0.9411021 18 13
#> 3: -8.123991 -0.9632126 -0.9718770 10 20
#> 4: -5.397180 -0.9192521 -0.9365610 11 11
#> 5: -9.785489 -0.9924813 -0.9949399 16 16
#> ---
#> 442: -2.502507 -0.8126749 -0.7662411 19 19
#> 443: -1.508060 -0.5524565 -0.5825105 20 20
#> 444: -5.061910 -0.9619721 -0.9606991 6 6
#> 445: -6.314742 -0.9658232 -0.9606853 17 17
#> 446: -17.018033 -0.9999538 -0.9999433 5 16
The scoring results has 10 columns:
-
varId: unique id of the variant as defined in the
id
column of the inputGRanges
orVRanges
object. If not defined, a custom unique identifier is generated in the format [seqname]:[range]:[ref_allele]>[alt_allele] - semId: the unique identifier of the SEM
- refSeq and altSeq: the optimally scored sequence for the reference and alternative alleles respectively
- refScore and altScore: the raw (unnormalized) binding affinity score for the reference and alternative alleles respectively
- refNorm and altNorm: the normalized binding affinity score for the reference and alternative alleles respectively
- refVarIndex and altVarIndex: The position of the frame’s starting indeces in the scored sequence
There are accessor functions to isolate each slot of the resulting data object:
# access the variants slot
variants(sempl_results)
#> VRanges object with 2 ranges and 5 metadata columns:
#> seqnames ranges strand ref alt totalDepth
#> <Rle> <IRanges> <Rle> <character> <characterOrRle> <integerOrRle>
#> [1] chr12 94136009 * G C <NA>
#> [2] chr19 10640062 * T A <NA>
#> refDepth altDepth sampleNames softFilterMatrix |
#> <integerOrRle> <integerOrRle> <factorOrRle> <matrix> |
#> [1] <NA> <NA> <NA> |
#> [2] <NA> <NA> <NA> |
#> id upstream downstream
#> <character> <character> <character>
#> [1] variant1 CTAATAAATACTTGTTGTT TTTAACAATAATGATGGTA
#> [2] variant2 GTAATCTTGGGCAAATTAT TAATCCTCTAAGGTCCACT
#> ref_seq alt_seq
#> <character> <character>
#> [1] CTAATAAATACTTGTTGTTG.. CTAATAAATACTTGTTGTTC..
#> [2] GTAATCTTGGGCAAATTATT.. GTAATCTTGGGCAAATTATA..
#> -------
#> seqinfo: 2 sequences from an unspecified genome; no seqlengths
#> hardFilters: NULL
# access the semData slot
semData(sempl_results)
#> Key: <SEM_KEY>
#> transcription_factor ensembl_id ebi_complex_ac uniprot_ac
#> <char> <char> <char> <char>
#> 1: TFAP2B ENSG00000008196 Q92481
#> 2: ARNT ENSG00000143437 P27540
#> 3: ATF1 ENSG00000123268 P18846
#> 4: ATF2 ENSG00000115966 P15336
#> 5: ATF3 ENSG00000162772 P18847
#> ---
#> 219: ZBTB7A ENSG00000178951 O95365
#> 220: ZFX ENSG00000005889 P17010
#> 221: ZNF281 ENSG00000162702 Q9Y2X9
#> 222: ZNF18 ENSG00000154957 P17022
#> 223: ZSCAN4 ENSG00000180532 Q8NAM6
#> PWM_id SEM SEM_baseline cell_type
#> <char> <char> <num> <char>
#> 1: AP2B_HUMAN.H11MO.0.B AP2B_HUMAN.SK-N-SH.sem -1.1610478 SK-N-SH
#> 2: ARNT_HUMAN.H11MO.0.B ARNT_HUMAN.GM12878.sem -1.8632604 GM12878
#> 3: ATF1_HUMAN.H11MO.0.B ATF1_HUMAN.K562.sem -2.9718855 K562
#> 4: ATF2_HUMAN.H11MO.0.B ATF2_HUMAN.HepG2.sem -1.4186930 HepG2
#> 5: ATF3_HUMAN.H11MO.0.A ATF3_HUMAN.HepG2.sem -2.1588687 HepG2
#> ---
#> 219: ZBT7A_HUMAN.H11MO.0.A ZBT7A_HUMAN.HepG2.sem -0.4056000 HepG2
#> 220: ZFX_HUMAN.H11MO.0.A ZFX_HUMAN.HepG2.sem -0.2478717 HepG2
#> 221: ZN281_HUMAN.H11MO.0.A ZN281_HUMAN.HepG2.sem -0.3926150 HepG2
#> 222: ZNF18_HUMAN.H11MO.0.C ZNF18_HUMAN.HEK293.sem -1.6459529 HEK293
#> 223: Zscan4_pwm_secondary ZSCAN4_secondary.sem -2.9110180 HEK293
#> neg_log10_pval chip_ENCODE_accession dnase_ENCODE_accession PWM_source
#> <num> <char> <char> <char>
#> 1: 6.476702 ENCFF303ORH ENCFF752OZB HOCOMOCOv11
#> 2: 12.081160 ENCFF538DNI ENCFF759OLD HOCOMOCOv11
#> 3: 56.056680 ENCFF885JPB ENCFF274YGF HOCOMOCOv11
#> 4: 27.821210 ENCFF899JAX ENCFF897NME HOCOMOCOv11
#> 5: 14.618530 ENCFF197JBM ENCFF897NME HOCOMOCOv11
#> ---
#> 219: 3.684286 ENCFF439DRT ENCFF897NME HOCOMOCOv11
#> 220: 3.443155 ENCFF329DCV ENCFF897NME HOCOMOCOv11
#> 221: 3.428189 ENCFF948PYK ENCFF897NME HOCOMOCOv11
#> 222: 3.793244 ENCFF888BHA ENCFF285OXK HOCOMOCOv11
#> 223: 11.973690 ENCFF485WHZ ENCFF910QHN UNIPROBE
#> SEM_KEY
#> <char>
#> 1: AP2B_HUMAN.SK-N-SH
#> 2: ARNT_HUMAN.GM12878
#> 3: ATF1_HUMAN.K562
#> 4: ATF2_HUMAN.HepG2
#> 5: ATF3_HUMAN.HepG2
#> ---
#> 219: ZBT7A_HUMAN.HepG2
#> 220: ZFX_HUMAN.HepG2
#> 221: ZN281_HUMAN.HepG2
#> 222: ZNF18_HUMAN.HEK293
#> 223: ZSCAN4_secondary
# access the scores slot
scores(sempl_results)
#> varId semId refSeq altSeq refScore
#> <char> <char> <char> <char> <num>
#> 1: variant1 AP2B_HUMAN.SK-N-SH GTTTAACAAT TCTTTAACAA -5.287583
#> 2: variant1 ARNT_HUMAN.GM12878 TTGTTTAAC TGTTGTTCT -5.048425
#> 3: variant1 ATF1_HUMAN.K562 ACTTGTTGTTG CTTTAACAATA -7.736531
#> 4: variant1 ATF2_HUMAN.HepG2 CTTGTTGTTGT CTTGTTGTTCT -5.049124
#> 5: variant1 ATF3_HUMAN.HepG2 TGTTGTTTAAC TGTTCTTTAAC -9.214165
#> ---
#> 442: variant2 ZBT7A_HUMAN.HepG2 TTTAATCCT TATAATCCT -2.821984
#> 443: variant2 ZFX_HUMAN.HepG2 TTAATCCTCT ATAATCCTCT -1.407772
#> 444: variant2 ZN281_HUMAN.HepG2 CTTGGGCAAATTATT CTTGGGCAAATTATA -5.109413
#> 445: variant2 ZNF18_HUMAN.HEK293 TATTTAATCCTC TATATAATCCTC -6.516791
#> 446: variant2 ZSCAN4_secondary TCTTGGGCAAATTATT TTATATAATCCTCTAA -17.314249
#> altScore refNorm altNorm refVarIndex altVarIndex
#> <num> <num> <num> <int> <int>
#> 1: -4.330368 -0.9427482 -0.8888423 20 19
#> 2: -5.948900 -0.8900564 -0.9411021 18 13
#> 3: -8.123991 -0.9632126 -0.9718770 10 20
#> 4: -5.397180 -0.9192521 -0.9365610 11 11
#> 5: -9.785489 -0.9924813 -0.9949399 16 16
#> ---
#> 442: -2.502507 -0.8126749 -0.7662411 19 19
#> 443: -1.508060 -0.5524565 -0.5825105 20 20
#> 444: -5.061910 -0.9619721 -0.9606991 6 6
#> 445: -6.314742 -0.9658232 -0.9606853 17 17
#> 446: -17.018033 -0.9999538 -0.9999433 5 16
We can view the optimal binding frame for any variant and SEM with
viewFrames
. Because alleles of a variant are scored
independently, the optimal binding site for each allele may be
different.
viewFrames(sempl_results, vid = "variant1", sid = "AP2B_HUMAN.SK-N-SH")
#> SEM: AP2B_HUMAN.SK-N-SH
#> Var: variant1
#> ref: CTAATAAATACTTGTTGTTGTTTAACAATAATGATGGTA
#> alt: CTAATAAATACTTGTTGTTCTTTAACAATAATGATGGTA
Visualization
plotSemMotifs(sempl_results,
variant = "variant2",
label = "transcription_factor")
plotSemVariants(sempl_results, semId = "TFAP4_HUMAN.HepG2")
Enrichment
# loading data
caqtl_file <- "caQTL_variants_overlappingPeaks_LD-r2-0.8_withLead.bed"
caqtls <- data.table::fread(file = caqtl_file, sep = "\t")[1:1000,]
proxy_ID_split <- strsplit(caqtls$proxy_ID, ":")
# formatting position and alleles
chr <- lapply(proxy_ID_split, `[[`, 1) |> unlist()
pos <- lapply(proxy_ID_split, `[[`, 2) |> unlist()
ref_allele <- lapply(proxy_ID_split, `[[`, 3) |> unlist()
alt_allele <- lapply(proxy_ID_split, `[[`, 4) |> unlist()
# building variant annotation obj
vr_e <- VariantAnnotation::VRanges(seqnames = chr,
ranges = pos,
ref = ref_allele,
alt = alt_allele)
data(sc)
# score variants
var_scores_e <- scoreVariants(vr = vr_e,
semList = sc,
bs_genome_obj=BSgenome.Hsapiens.UCSC.hg38::Hsapiens)
plotSemVariants(var_scores_e, semId = "MA0151.1")
plotSemMotifs(var_scores_e, "chr1:906982:C>T", label = "transcription_factor")
Binding Enrichment
Given this is a large set of variants, we might want to check if some transcription factors are enriched for binding in one allele or the other.
Here, we will test if any of these transcription factors are bound more than expected in the reference alleles of this set.
gr_e <- GenomicRanges::GRanges(seqnames = chr,
ranges = pos,
allele = ref_allele)
e <- enrichSEMs(x = gr_e, semList = sc,
bs_genome_obj = BSgenome.Hsapiens.UCSC.hg38::Hsapiens)
plt <- plotEnrich(e, semList = sc)
plt
Extras
Loading SEMs
First, we create a list of file paths to the .sem files we want to
include. We will also load the meta data for each sem in a
data.table
object. If meta data is used, all SEMs must be
represented in the meta data table.
sem_files <- list.files("SEMs/",
pattern = ".sem", full.names = TRUE)
metadata_file <- "SEMs/sempl_metadata.csv"
meta <- data.table::fread(file = metadata_file, sep = ",")
head(meta)
We will load the matrix and meta data for all SEMs in a single
SNPEffectMatrixCollection
object. This object has two
slots, one containing a named list of the matrices and a second slot
containing our meta data table with a key column connecting the meta
data to the names of the matrices.
sc <- loadSEMCollection(semFiles = sem_files,
semMetaData = meta,
semMetaKey = "SEM")
sc
devtools::session_info()
#> ─ Session info ───────────────────────────────────────────────────────────────
#> setting value
#> version R version 4.4.1 (2024-06-14)
#> os macOS Monterey 12.6
#> system aarch64, darwin20
#> ui X11
#> language en
#> collate en_US.UTF-8
#> ctype en_US.UTF-8
#> tz America/New_York
#> date 2025-07-30
#> pandoc 3.1.1 @ /Applications/RStudio.app/Contents/Resources/app/quarto/bin/tools/ (via rmarkdown)
#>
#> ─ Packages ───────────────────────────────────────────────────────────────────
#> package * version date (UTC) lib source
#> abind 1.4-5 2016-07-21 [1] CRAN (R 4.4.0)
#> AnnotationDbi 1.66.0 2024-06-16 [1] Bioconductor 3.19 (R 4.4.0)
#> Biobase * 2.64.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> BiocGenerics * 0.50.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> BiocIO 1.14.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> BiocParallel 1.38.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> Biostrings * 2.72.1 2024-06-02 [1] Bioconductor 3.19 (R 4.4.0)
#> bit 4.0.5 2022-11-15 [1] CRAN (R 4.4.0)
#> bit64 4.0.5 2020-08-30 [1] CRAN (R 4.4.0)
#> bitops 1.0-8 2024-07-29 [1] CRAN (R 4.4.0)
#> blob 1.2.4 2023-03-17 [1] CRAN (R 4.4.0)
#> BSgenome 1.72.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> BSgenome.Hsapiens.UCSC.hg38 1.4.5 2025-01-07 [1] Bioconductor
#> bslib 0.8.0 2024-07-29 [1] CRAN (R 4.4.0)
#> cachem 1.1.0 2024-05-16 [1] CRAN (R 4.4.0)
#> cli 3.6.3 2024-06-21 [1] CRAN (R 4.4.0)
#> codetools 0.2-20 2024-03-31 [1] CRAN (R 4.4.1)
#> colorspace 2.1-1 2024-07-26 [1] CRAN (R 4.4.0)
#> crayon 1.5.3 2024-06-20 [1] CRAN (R 4.4.0)
#> curl 5.2.1 2024-03-01 [1] CRAN (R 4.4.0)
#> data.table 1.15.4 2024-03-30 [1] CRAN (R 4.4.0)
#> DBI 1.2.3 2024-06-02 [1] CRAN (R 4.4.0)
#> DelayedArray 0.30.1 2024-05-30 [1] Bioconductor 3.19 (R 4.4.0)
#> desc 1.4.3 2023-12-10 [1] CRAN (R 4.4.0)
#> devtools 2.4.5 2022-10-11 [1] CRAN (R 4.4.0)
#> digest 0.6.36 2024-06-23 [1] CRAN (R 4.4.0)
#> dplyr 1.1.4 2023-11-17 [1] CRAN (R 4.4.0)
#> ellipsis 0.3.2 2021-04-29 [1] CRAN (R 4.4.0)
#> evaluate 0.24.0 2024-06-10 [1] CRAN (R 4.4.0)
#> fansi 1.0.6 2023-12-08 [1] CRAN (R 4.4.0)
#> farver 2.1.2 2024-05-13 [1] CRAN (R 4.4.0)
#> fastmap 1.2.0 2024-05-15 [1] CRAN (R 4.4.0)
#> fs 1.6.4 2024-04-25 [1] CRAN (R 4.4.0)
#> generics 0.1.3 2022-07-05 [1] CRAN (R 4.4.0)
#> GenomeInfoDb * 1.40.1 2024-06-16 [1] Bioconductor 3.19 (R 4.4.0)
#> GenomeInfoDbData 1.2.12 2024-08-21 [1] Bioconductor
#> GenomicAlignments 1.40.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> GenomicFeatures 1.56.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> GenomicRanges * 1.56.1 2024-06-16 [1] Bioconductor 3.19 (R 4.4.0)
#> ggplot2 3.5.1 2024-04-23 [1] CRAN (R 4.4.0)
#> ggrepel 0.9.5 2024-01-10 [1] CRAN (R 4.4.0)
#> glue 1.7.0 2024-01-09 [1] CRAN (R 4.4.0)
#> gtable 0.3.5 2024-04-22 [1] CRAN (R 4.4.0)
#> highr 0.11 2024-05-26 [1] CRAN (R 4.4.0)
#> htmltools 0.5.8.1 2024-04-04 [1] CRAN (R 4.4.0)
#> htmlwidgets 1.6.4 2023-12-06 [1] CRAN (R 4.4.0)
#> httpuv 1.6.15 2024-03-26 [1] CRAN (R 4.4.0)
#> httr 1.4.7 2023-08-15 [1] CRAN (R 4.4.0)
#> IRanges * 2.38.1 2024-07-03 [1] Bioconductor 3.19 (R 4.4.1)
#> jquerylib 0.1.4 2021-04-26 [1] CRAN (R 4.4.0)
#> jsonlite 1.8.8 2023-12-04 [1] CRAN (R 4.4.0)
#> KEGGREST 1.44.1 2024-06-19 [1] Bioconductor 3.19 (R 4.4.0)
#> knitr 1.48 2024-07-07 [1] CRAN (R 4.4.0)
#> labeling 0.4.3 2023-08-29 [1] CRAN (R 4.4.0)
#> later 1.3.2 2023-12-06 [1] CRAN (R 4.4.0)
#> lattice 0.22-6 2024-03-20 [1] CRAN (R 4.4.1)
#> lifecycle 1.0.4 2023-11-07 [1] CRAN (R 4.4.0)
#> magrittr 2.0.3 2022-03-30 [1] CRAN (R 4.4.0)
#> Matrix 1.7-0 2024-04-26 [1] CRAN (R 4.4.1)
#> MatrixGenerics * 1.16.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> matrixStats * 1.3.0 2024-04-11 [1] CRAN (R 4.4.0)
#> memoise 2.0.1 2021-11-26 [1] CRAN (R 4.4.0)
#> mime 0.12 2021-09-28 [1] CRAN (R 4.4.0)
#> miniUI 0.1.1.1 2018-05-18 [1] CRAN (R 4.4.0)
#> munsell 0.5.1 2024-04-01 [1] CRAN (R 4.4.0)
#> pillar 1.9.0 2023-03-22 [1] CRAN (R 4.4.0)
#> pkgbuild 1.4.4 2024-03-17 [1] CRAN (R 4.4.0)
#> pkgconfig 2.0.3 2019-09-22 [1] CRAN (R 4.4.0)
#> pkgdown 2.1.0 2024-07-06 [1] CRAN (R 4.4.0)
#> pkgload 1.4.0 2024-06-28 [1] CRAN (R 4.4.0)
#> png 0.1-8 2022-11-29 [1] CRAN (R 4.4.0)
#> profvis 0.3.8 2023-05-02 [1] CRAN (R 4.4.0)
#> promises 1.3.0 2024-04-05 [1] CRAN (R 4.4.0)
#> purrr 1.0.2 2023-08-10 [1] CRAN (R 4.4.0)
#> R6 2.5.1 2021-08-19 [1] CRAN (R 4.4.0)
#> ragg 1.3.2 2024-05-15 [1] CRAN (R 4.4.0)
#> Rcpp 1.0.13 2024-07-17 [1] CRAN (R 4.4.0)
#> RCurl 1.98-1.16 2024-07-11 [1] CRAN (R 4.4.0)
#> remotes 2.5.0 2024-03-17 [1] CRAN (R 4.4.0)
#> restfulr 0.0.15 2022-06-16 [1] CRAN (R 4.4.0)
#> rjson 0.2.22 2024-08-20 [1] CRAN (R 4.4.1)
#> rlang 1.1.4 2024-06-04 [1] CRAN (R 4.4.0)
#> rmarkdown 2.28 2024-08-17 [1] CRAN (R 4.4.0)
#> Rsamtools * 2.20.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> RSQLite 2.3.7 2024-05-27 [1] CRAN (R 4.4.0)
#> rstudioapi 0.16.0 2024-03-24 [1] CRAN (R 4.4.0)
#> rtracklayer 1.64.0 2024-05-06 [1] Bioconductor 3.19 (R 4.4.0)
#> S4Arrays 1.4.1 2024-05-30 [1] Bioconductor 3.19 (R 4.4.0)
#> S4Vectors * 0.42.1 2024-07-03 [1] Bioconductor 3.19 (R 4.4.1)
#> sass 0.4.9 2024-03-15 [1] CRAN (R 4.4.0)
#> scales 1.3.0 2023-11-28 [1] CRAN (R 4.4.0)
#> SEMplR * 0.0.0.9000 2025-06-24 [1] local
#> sessioninfo 1.2.2 2021-12-06 [1] CRAN (R 4.4.0)
#> shiny 1.9.1 2024-08-01 [1] CRAN (R 4.4.0)
#> SparseArray 1.4.8 2024-05-30 [1] Bioconductor 3.19 (R 4.4.0)
#> stringi 1.8.4 2024-05-06 [1] CRAN (R 4.4.0)
#> stringr 1.5.1 2023-11-14 [1] CRAN (R 4.4.0)
#> SummarizedExperiment * 1.34.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> systemfonts 1.1.0 2024-05-15 [1] CRAN (R 4.4.0)
#> textshaping 0.4.0 2024-05-24 [1] CRAN (R 4.4.0)
#> tibble 3.2.1 2023-03-20 [1] CRAN (R 4.4.0)
#> tidyselect 1.2.1 2024-03-11 [1] CRAN (R 4.4.0)
#> UCSC.utils 1.0.0 2024-05-06 [1] Bioconductor 3.19 (R 4.4.0)
#> urlchecker 1.0.1 2021-11-30 [1] CRAN (R 4.4.0)
#> usethis 3.0.0 2024-07-29 [1] CRAN (R 4.4.0)
#> utf8 1.2.4 2023-10-22 [1] CRAN (R 4.4.0)
#> VariantAnnotation * 1.48.1 2023-11-15 [1] Bioconductor
#> vctrs 0.6.5 2023-12-01 [1] CRAN (R 4.4.0)
#> withr 3.0.1 2024-07-31 [1] CRAN (R 4.4.0)
#> xfun 0.49 2024-10-31 [1] CRAN (R 4.4.1)
#> XML 3.99-0.17 2024-06-25 [1] CRAN (R 4.4.0)
#> xtable 1.8-4 2019-04-21 [1] CRAN (R 4.4.0)
#> XVector * 0.44.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#> yaml 2.3.10 2024-07-26 [1] CRAN (R 4.4.0)
#> zlibbioc 1.50.0 2024-04-30 [1] Bioconductor 3.19 (R 4.4.0)
#>
#> [1] /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/library
#>
#> ──────────────────────────────────────────────────────────────────────────────