EvoGeneX is a package to infer mode of gene expression evolution from expression data.
library(tidyverse)
library(EvoGeneX)
library(png)
library(knitr)
library(ape)
In addition to the expression data, EvoGeneX needs 1) a phylogenetic tree, and 2) information about different regimes operating on different parts of the phylogenetic tree.
The phylogenetic tree should have exactly the set of species for which gene expression data is collected as the leaf nodes. It must have numeric length on the branches abstracting evolutionary time since the root denoting the most recent ancestor of the species under study. It should be specified in Newick format. A copy of the nine-species Drosphila phylogeny used in the EvoGeneX paper is kept at inst/extdata/drosophila9.newick
.
dros9_newick_file <- system.file("extdata", "drosophila9.newick", package = "EvoGeneX")
phylo_tree <- read.tree(dros9_newick_file)
plot(phylo_tree)
If EvoGeneX is used to infer about adaptation of expression towards different parts of the phylogenetic tree, it must be specified as a data frame containing three columns each taking string values: 1) node
, 2) node2
and 3) regime
and as many rows as the nodes in the tree including the non-leaf internal nodes. For the leaf nodes, node
gives species name, regime
gives a name for the regime and node2
is kept empty. Each internal node x
is specified by two leaf species node
and node2
such that x
is the most recent common ancestor of node
and node2
. Note that x
can be specified using different combinations of node
and node2
, use any one of them.
If EvoGeneX is used to infer if a gene evolved in a constrained way, not necessarily in an adaptive way, all nodes in the tree are assumed to be in a single regime. Nevertheless, and we need to specify that in a file as given in inst/extdata/regime_global.csv
.
single_regime_file <- system.file("extdata", "regime_global.csv", package = "EvoGeneX")
single_regime <- read.csv(single_regime_file)
single_regime
## node node2 regime
## 1 dmel global
## 2 dyak global
## 3 dana global
## 4 dpse global
## 5 dper global
## 6 dwil global
## 7 dvir global
## 8 dmoj global
## 9 dgri global
## 10 dmel dyak global
## 11 dmel dana global
## 12 dpse dper global
## 13 dmel dper global
## 14 dmel dwil global
## 15 dvir dmoj global
## 16 dvir dgri global
## 17 dmel dgri global
EvoGeneX needs exactly nspecies
x nreplicate
expression values where nspecies
and nreplicate
denote the number of species and replicates, respectively. Currently, EvoGeneX does not support missing values of a few replicates of some of the species, which may be supported in a future version. These nspecies
x nreplicate
expression values can be provided as a R data.frame
either in the matrix (wide
) or vector format (tall
).
In the wide
format the data.frame
must contain a column for the species and the rest of the columns will be treated as replicates. The default name of the species column is species
though it can be changed through the species_col
parameter of the fit
function (see below). There must be exactly nspecies
number of rows. An example expression data for 9 Drosophila species in 4 replicates is given in inst/extdata/HD_M_FBgn0000008.csv
.
sample_data_file <- system.file("extdata", "HD_M_FBgn0000008.csv", package = "EvoGeneX")
wide <- read.csv(sample_data_file)
wide
## species R1 R2 R3 R4
## 1 dana 204.5885 148.66722 230.49887 199.3192
## 2 dgri 318.4320 258.96714 240.11517 259.6903
## 3 dmel 241.2408 212.74946 204.11434 190.5641
## 4 dmoj 132.7657 69.73923 65.77837 110.6977
## 5 dper 107.1932 109.08483 124.48701 140.8748
## 6 dpse 245.1403 188.31841 228.28671 211.3906
## 7 dvir 543.9844 696.80922 663.95690 785.1595
## 8 dwil 155.7480 214.93546 282.80718 316.0916
## 9 dyak 106.0654 117.31824 102.06957 108.0820
In the tall
format the data.frame
must contain at least 3 columns for the species, replicates and expression values and the rest of the columns will simply be ignored. The default name of the 3 columns are species
, replicate
and exprval
, respectively, though they can be changed through the parameters species_col
, replicate_col
and exprval_col
, respectively, of the fit
function. Note that the parameters replicate_col
and exprval_col
are not relevant in wide
format. There must be exactly nspecies
x nreplicate
number of rows in tall
format. Following code shows the same example expression data in tall
format.
tall <- wide %>% gather("replicate", "exprval", -species)
tall
## species replicate exprval
## 1 dana R1 204.58850
## 2 dgri R1 318.43200
## 3 dmel R1 241.24080
## 4 dmoj R1 132.76570
## 5 dper R1 107.19320
## 6 dpse R1 245.14030
## 7 dvir R1 543.98440
## 8 dwil R1 155.74800
## 9 dyak R1 106.06540
## 10 dana R2 148.66722
## 11 dgri R2 258.96714
## 12 dmel R2 212.74946
## 13 dmoj R2 69.73923
## 14 dper R2 109.08483
## 15 dpse R2 188.31841
## 16 dvir R2 696.80922
## 17 dwil R2 214.93546
## 18 dyak R2 117.31824
## 19 dana R3 230.49887
## 20 dgri R3 240.11517
## 21 dmel R3 204.11434
## 22 dmoj R3 65.77837
## 23 dper R3 124.48701
## 24 dpse R3 228.28671
## 25 dvir R3 663.95690
## 26 dwil R3 282.80718
## 27 dyak R3 102.06957
## 28 dana R4 199.31920
## 29 dgri R4 259.69030
## 30 dmel R4 190.56410
## 31 dmoj R4 110.69770
## 32 dper R4 140.87480
## 33 dpse R4 211.39060
## 34 dvir R4 785.15950
## 35 dwil R4 316.09160
## 36 dyak R4 108.08200
Now we show how to infer whether the expression of the gene in the example data evolved neutrally or in a constrained manner. The example script examples/example_constrained.R
gives the complete code used in this section.
evog <- EvoGeneX()
evog$setTree(dros9_newick_file)
evog$setRegimes(single_regime_file)
ou_res <- evog$fit(wide, alpha = 0.1, gamma = 0.01)
print(ou_res)
## $optim.diagn
## $optim.diagn$convergence
## [1] 1
##
## $optim.diagn$message
## NULL
##
##
## $theta
## global
## 247.6608
##
## $alpha
## [1] 13.76675
##
## $sigma.sq
## [1] 785188.3
##
## $gamma
## [1] 0.002799634
##
## $loglik
## [1] -206.8637
We decide whether a gene has evolved neutrally or in a constrained manner by fitting the data to both replicated models, Ornstein-Uhlenbeck and Brownian Motion and checking which model is most likely supported by the data. So we fit the data to the replicated Brownian Motion model as well.
brown <- Brown()
brown$setTree(dros9_newick_file)
brown_res <- brown$fit(wide, gamma = 0.01)
print(brown_res)
## $optim.diagn
## $optim.diagn$convergence
## [1] 1
##
## $optim.diagn$message
## NULL
##
##
## $theta
## [1] 285.5768
##
## $sigma.sq
## [1] 167801.1
##
## $gamma
## [1] 0.01346066
##
## $loglik
## [1] -209.7939
We use log-likelihood ratio test to check which model most likely support the data. That is implemented by \(\chi^2\) test on the test statistics -2 * (log-likelihood in null model - loglikelihood in alternate model)) taking Brownian motion as null and OU as alternate.
# degrees of freedom under different models
ou_dof <- (
1 # alpha
+ 1 # sigma.sq
+ 1 # theta
+ 1 # gamma
)
brown_dof <- (
1 # sigma.sq
+ 1 # theta
+ 1 # gamma
)
# loglikelihood ratio test OU VS Brownian motion
pvalue <- 1 - pchisq((ou_res$loglik - brown_res$loglik) * 2,
(ou_dof - brown_dof))
constrained_vs_neutral <- ifelse(pvalue < 0.05, "constrained", "neutral")
constrained_vs_neutral
## [1] "constrained"
Now we show how to infer whether the expression of the gene in the example data evolved in an adaptive manner. The example script examples/example_adaptive.R
gives the complete code used in this section.
To test if a gene’s expression evolved in an adaptive manner, we need to specify at least two regimes and which of the nodes in the phylogenetic tree including the non-leaf internal nodes belong to each of the regimes. Here we consider a two-regime example where the species in Drosophila and Sophophora sub-genra belong to different regimes.
two_regime_image <- system.file("extdata", "phylogeny_two_regime.png", package = "EvoGeneX")
knitr::include_graphics(two_regime_image)
## Warning in knitr::include_graphics(two_regime_image): It is highly recommended
## to use relative paths for images. You had absolute paths: "/Users/pals2/Library/
## R/x86_64/4.1/library/EvoGeneX/extdata/phylogeny_two_regime.png"
two_regime_file <- system.file("extdata", "regime_tworegime.csv", package = "EvoGeneX")
two_regime <- read.csv(two_regime_file)
two_regime
## node node2 regime
## 1 dmel sophophora
## 2 dyak sophophora
## 3 dana sophophora
## 4 dpse sophophora
## 5 dper sophophora
## 6 dwil sophophora
## 7 dvir drosophila
## 8 dmoj drosophila
## 9 dgri drosophila
## 10 dmel dyak sophophora
## 11 dmel dana sophophora
## 12 dpse dper sophophora
## 13 dmel dper sophophora
## 14 dmel dwil sophophora
## 15 dvir dmoj drosophila
## 16 dvir dgri drosophila
## 17 dmel dgri sophophora
Like constrained case, we call a gene adaptive if both the neutral Bownian Motion model and constrained model (OU with a single optimum expression) are rejected in favor of adaptive evolution model (OU with multiple optima).
evog2 <- EvoGeneX()
evog2$setTree(dros9_newick_file)
evog2$setRegimes(two_regime_file)
ou2_res <- evog$fit(wide, alpha = 0.1, gamma = 0.01)
print(ou2_res)
## $optim.diagn
## $optim.diagn$convergence
## [1] 1
##
## $optim.diagn$message
## NULL
##
##
## $theta
## global
## 247.6608
##
## $alpha
## [1] 13.76675
##
## $sigma.sq
## [1] 785188.3
##
## $gamma
## [1] 0.002799634
##
## $loglik
## [1] -206.8637
# degrees of freedom under different models
ou2_dof <- (
1 # alpha
+ 1 # sigma.sq
+ 2 # theta
+ 1 # gamma
)
# loglikelihood ratio test OU2 VS replicated Brownian Motion
ou2_vs_bm_pvalue <- 1 - pchisq((ou2_res$loglik - brown_res$loglik) * 2,
(ou2_dof - brown_dof))
# loglikelihood ratio test OU2 VS OU
ou2_vs_ou_pvalue <- 1 - pchisq((ou2_res$loglik - ou_res$loglik) * 2,
(ou2_dof - ou_dof))
# gene is adaptive if both OU2vsBM and OU2vsOU tests are successfull
is_adaptive <- ifelse(max(ou2_vs_bm_pvalue, ou2_vs_ou_pvalue) < 0.05,
"adaptive", "not-adaptive")
is_adaptive
## [1] "not-adaptive"
In the main paper for EvoGeneX, we apply EvoGeneX to replicated gene expression (bulk RNA-seq) data from 5 different body-parts of both sexes from 9 Drosophila species. The EvoGeneX package includes this Drosophila data which can be accessed using function.
# Describe the data in the package
data(package="EvoGeneX")
# Access one example body-part sample
data(FemaleHead, package = "EvoGeneX")
head(FemaleHead)
## FBgnDmel SymbolDmel F_HD_dmel_R1 F_HD_dmel_R2 F_HD_dmel_R3 F_HD_dmel_R4
## 1 FBgn0000008 a 213.7328471 195.40689 203.6376820 224.4284
## 2 FBgn0000014 abd-A 0.0000000 0.00000 0.8970823 0.0000
## 3 FBgn0000015 Abd-B 0.6399187 0.00000 0.4485411 0.0000
## 4 FBgn0000017 Abl 1114.6018054 1023.84673 1020.2657126 867.4850
## 5 FBgn0000018 abo 72.3108136 65.82127 61.8986787 68.1024
## 6 FBgn0000022 ac 0.0000000 0.00000 0.0000000 0.0000
## F_HD_dyak_R1 F_HD_dyak_R2 F_HD_dyak_R3 F_HD_dyak_R4 F_HD_dana_R1 F_HD_dana_R2
## 1 83.264808 83.28962 106.178657 86.97742 192.11953 182.5986
## 2 0.000000 0.00000 0.000000 0.00000 0.00000 0.0000
## 3 1.460786 0.00000 1.106028 0.00000 0.00000 0.0000
## 4 679.265541 672.09424 569.972931 430.40372 516.16113 486.0537
## 5 45.284369 35.62677 30.968775 17.03681 38.42391 60.4283
## 6 0.000000 0.00000 0.000000 0.00000 0.00000 0.0000
## F_HD_dana_R3 F_HD_dana_R4 F_HD_dpse_R1 F_HD_dpse_R2 F_HD_dpse_R3 F_HD_dpse_R4
## 1 232.88030 179.3616 254.88975 230.737764 224.5431265 239.461144
## 2 0.00000 0.0000 0.00000 2.996594 0.0000000 0.000000
## 3 0.00000 0.0000 0.00000 0.000000 0.7906448 3.005788
## 4 493.74235 390.6436 930.23867 1018.842076 952.7269980 891.717230
## 5 40.61866 54.7205 54.46362 37.956862 47.4386887 49.094544
## 6 0.00000 0.0000 0.00000 0.000000 0.0000000 0.000000
## F_HD_dper_R1 F_HD_dper_R2 F_HD_dper_R3 F_HD_dper_R4 F_HD_dwil_R1 F_HD_dwil_R2
## 1 323.30107 410.05403 240.198864 173.9355335 273.9932826 262.28710
## 2 0.00000 0.00000 2.151547 0.5532203 0.9546804 0.00000
## 3 0.00000 0.00000 0.000000 0.0000000 0.0000000 0.00000
## 4 437.08406 495.45008 464.016870 401.0846832 1051.1031505 904.11905
## 5 34.68921 36.88041 40.162202 48.1301620 46.7793409 43.20023
## 6 0.00000 0.00000 0.000000 0.0000000 0.0000000 0.00000
## F_HD_dwil_R3 F_HD_dwil_R4 F_HD_dvir_R1 F_HD_dvir_R2 F_HD_dvir_R3 F_HD_dvir_R4
## 1 360.85678 293.994596 661.6476906 723.3232 634.99492 577.25510
## 2 0.00000 1.185462 0.7031325 0.0000 0.00000 0.00000
## 3 0.00000 0.000000 0.0000000 0.0000 0.00000 0.00000
## 4 1153.82648 1038.464784 544.2245617 607.5041 641.47446 527.85038
## 5 40.53102 50.974870 87.1884311 81.9475 73.70477 79.74094
## 6 0.00000 0.000000 0.0000000 0.0000 0.00000 0.00000
## F_HD_dmoj_R1 F_HD_dmoj_R2 F_HD_dmoj_R3 F_HD_dmoj_R4 F_HD_dgri_R1 F_HD_dgri_R2
## 1 74.647080 85.21123 87.1195640 71.14443 297.510683 200.258298
## 2 0.000000 0.00000 0.7992621 0.00000 0.000000 0.000000
## 3 0.829412 0.00000 0.7992621 0.00000 1.662071 1.982755
## 4 624.547237 675.83562 557.0856525 558.99191 741.283601 695.947155
## 5 112.800032 115.13273 97.5099707 91.47140 63.158692 69.396440
## 6 0.000000 0.00000 0.0000000 0.00000 0.000000 0.000000
## F_HD_dgri_R3 F_HD_dgri_R4
## 1 274.67882 248.98749
## 2 0.00000 0.00000
## 3 0.00000 0.00000
## 4 858.97368 594.83646
## 5 65.05551 61.20695
## 6 0.00000 0.00000
Here we show how we can check multiple genes for constrained expression evolution using Drosophila dataset. Similar can be done for adaptive evolution.
# take 5 genes from FemaleHead dataset
df <- FemaleHead[1:5,]
# setup EvoGeneX object
evog <- EvoGeneX()
evog$setTree(dros9_newick_file)
evog$setRegimes(single_regime_file)
# setup Brown object
brown <- Brown()
brown$setTree(dros9_newick_file)
# degrees of freedom under different models
ou_dof <- (
1 # alpha
+ 1 # sigma.sq
+ 1 # theta
+ 1 # gamma
)
brown_dof <- (
1 # sigma.sq
+ 1 # theta
+ 1 # gamma
)
fdr_cutoff <- 0.05
process_single_gene <- function(data_wide) {
# process single gene
# all 36 expression values are in a single row
# need convert it to either 'wide' or 'tall' format
# we convert into 'tall' format
data_tall <- (
data_wide
%>% gather("sample", "exprval")
%>% separate(sample, c("sex", "tissue", "species", "replicate"))
%>% select("species", "replicate", "exprval")
)
ou_res <- evog$fit(data_tall, format = "tall", alpha = 0.1, gamma = 0.01)
brown_res <- brown$fit(data_tall, format = "tall", gamma = 0.01)
# loglikelihood ratio test EvoGeneX VS replicated Brownian motion
pvalue <- 1 - pchisq((ou_res$loglik - brown_res$loglik) * 2,
(ou_dof - brown_dof))
}
res <- (
df
%>% group_by(FBgnDmel, SymbolDmel)
%>% summarize(pvalue = process_single_gene(cur_data()), .groups = "keep")
%>% ungroup()
%>% mutate(qvalue = p.adjust(pvalue, method = "fdr"))
%>% mutate(constrained_vs_neutral = ifelse(qvalue < fdr_cutoff,
"constrained", "neutral"))
)
citation("EvoGeneX")
##
## To cite EvoGeneX in publications use:
##
## Soumitra Pal, Brian C. Oliver, and Teresa M. Przytycka. Modeling gene
## expression evolution with EvoGeneX uncovers differences in evolution
## of species, organs and sexes. bioRxiv (2021): 2020-01.
## <10.1101/2020.01.06.895615>.
##
## A BibTeX entry for LaTeX users is
##
## @Article{,
## title = {Stochastic modeling of gene expression evolution uncovers tissue and sex specific properties of expression evolution in the Drosophila genus},
## author = {Soumitra Pal and Brian C. Oliver and Teresa M. Przytycka},
## journal = {Journal of Computational Biology},
## year = {2022},
## volume = {xx},
## number = {yy},
## pages = {ppp},
## doi = {10.1101/2020.01.06.895615},
## }
devtools::session_info()
## ─ Session info ───────────────────────────────────────────────────────────────
## setting value
## version R version 4.1.2 (2021-11-01)
## os macOS Catalina 10.15.7
## system x86_64, darwin17.0
## ui X11
## language (EN)
## collate en_US.UTF-8
## ctype en_US.UTF-8
## tz America/New_York
## date 2022-06-03
## pandoc 2.3.1 @ /Applications/RStudio.app/Contents/MacOS/pandoc/ (via rmarkdown)
##
## ─ Packages ───────────────────────────────────────────────────────────────────
## package * version date (UTC) lib source
## ape * 5.6-2 2022-03-02 [1] CRAN (R 4.1.2)
## assertthat 0.2.1 2019-03-21 [1] CRAN (R 4.1.0)
## backports 1.4.1 2021-12-13 [1] CRAN (R 4.1.0)
## brio 1.1.3 2021-11-30 [1] CRAN (R 4.1.0)
## broom 0.7.12 2022-01-28 [1] CRAN (R 4.1.2)
## bslib 0.3.1 2021-10-06 [1] CRAN (R 4.1.0)
## cachem 1.0.6 2021-08-19 [1] CRAN (R 4.1.0)
## callr 3.7.0 2021-04-20 [1] CRAN (R 4.1.0)
## cellranger 1.1.0 2016-07-27 [1] CRAN (R 4.1.0)
## cli 3.2.0 2022-02-14 [1] CRAN (R 4.1.2)
## colorspace 2.0-3 2022-02-21 [1] CRAN (R 4.1.2)
## crayon 1.5.1 2022-03-26 [1] CRAN (R 4.1.2)
## DBI 1.1.2 2021-12-20 [1] CRAN (R 4.1.0)
## dbplyr 2.1.1 2021-04-06 [1] CRAN (R 4.1.0)
## desc 1.4.1 2022-03-06 [1] CRAN (R 4.1.2)
## devtools 2.4.3 2021-11-30 [1] CRAN (R 4.1.0)
## digest 0.6.29 2021-12-01 [1] CRAN (R 4.1.0)
## dplyr * 1.0.8 2022-02-08 [1] CRAN (R 4.1.2)
## ellipsis 0.3.2 2021-04-29 [1] CRAN (R 4.1.0)
## evaluate 0.15 2022-02-18 [1] CRAN (R 4.1.2)
## EvoGeneX * 0.9.9.0 2022-06-01 [1] Github (ncbi/EvoGeneX@f7db424)
## fansi 1.0.3 2022-03-24 [1] CRAN (R 4.1.2)
## fastmap 1.1.0 2021-01-25 [1] CRAN (R 4.1.0)
## forcats * 0.5.1 2021-01-27 [1] CRAN (R 4.1.0)
## fs 1.5.2 2021-12-08 [1] CRAN (R 4.1.0)
## generics 0.1.2 2022-01-31 [1] CRAN (R 4.1.2)
## ggplot2 * 3.3.5 2021-06-25 [1] CRAN (R 4.1.0)
## glue 1.6.2 2022-02-24 [1] CRAN (R 4.1.2)
## gtable 0.3.0 2019-03-25 [1] CRAN (R 4.1.0)
## haven 2.4.3 2021-08-04 [1] CRAN (R 4.1.0)
## highr 0.9 2021-04-16 [1] CRAN (R 4.1.0)
## hms 1.1.1 2021-09-26 [1] CRAN (R 4.1.0)
## htmltools 0.5.2 2021-08-25 [1] CRAN (R 4.1.0)
## httr 1.4.2 2020-07-20 [1] CRAN (R 4.1.0)
## jquerylib 0.1.4 2021-04-26 [1] CRAN (R 4.1.0)
## jsonlite 1.8.0 2022-02-22 [1] CRAN (R 4.1.2)
## knitr * 1.38 2022-03-25 [1] CRAN (R 4.1.2)
## lattice 0.20-45 2021-09-22 [2] CRAN (R 4.1.2)
## lifecycle 1.0.1 2021-09-24 [1] CRAN (R 4.1.0)
## lubridate 1.8.0 2021-10-07 [1] CRAN (R 4.1.0)
## magrittr 2.0.3 2022-03-30 [1] CRAN (R 4.1.2)
## Matrix 1.3-4 2021-06-01 [2] CRAN (R 4.1.2)
## memoise 2.0.1 2021-11-26 [1] CRAN (R 4.1.0)
## modelr 0.1.8 2020-05-19 [1] CRAN (R 4.1.0)
## munsell 0.5.0 2018-06-12 [1] CRAN (R 4.1.0)
## nlme 3.1-153 2021-09-07 [2] CRAN (R 4.1.2)
## nloptr * 2.0.0 2022-01-26 [1] CRAN (R 4.1.2)
## ouch * 2.17 2021-05-16 [1] CRAN (R 4.1.0)
## pillar 1.7.0 2022-02-01 [1] CRAN (R 4.1.2)
## pkgbuild 1.3.1 2021-12-20 [1] CRAN (R 4.1.0)
## pkgconfig 2.0.3 2019-09-22 [1] CRAN (R 4.1.0)
## pkgload 1.2.4 2021-11-30 [1] CRAN (R 4.1.0)
## plyr 1.8.7 2022-03-24 [1] CRAN (R 4.1.2)
## png * 0.1-7 2013-12-03 [1] CRAN (R 4.1.0)
## prettyunits 1.1.1 2020-01-24 [1] CRAN (R 4.1.0)
## processx 3.5.3 2022-03-25 [1] CRAN (R 4.1.2)
## ps 1.6.0 2021-02-28 [1] CRAN (R 4.1.0)
## purrr * 0.3.4 2020-04-17 [1] CRAN (R 4.1.0)
## R6 2.5.1 2021-08-19 [1] CRAN (R 4.1.0)
## Rcpp * 1.0.8.3 2022-03-17 [1] CRAN (R 4.1.2)
## RcppEigen * 0.3.3.9.1 2020-12-17 [1] CRAN (R 4.1.0)
## readr * 2.1.2 2022-01-30 [1] CRAN (R 4.1.2)
## readxl 1.4.0 2022-03-28 [1] CRAN (R 4.1.2)
## remotes 2.4.2 2021-11-30 [1] CRAN (R 4.1.0)
## reprex 2.0.1 2021-08-05 [1] CRAN (R 4.1.0)
## reshape2 1.4.4 2020-04-09 [1] CRAN (R 4.1.0)
## rlang 1.0.2 2022-03-04 [1] CRAN (R 4.1.2)
## rmarkdown 2.13 2022-03-10 [1] CRAN (R 4.1.2)
## rprojroot 2.0.3 2022-04-02 [1] CRAN (R 4.1.2)
## rstudioapi 0.13 2020-11-12 [1] CRAN (R 4.1.0)
## rvest 1.0.2 2021-10-16 [1] CRAN (R 4.1.0)
## sass 0.4.1 2022-03-23 [1] CRAN (R 4.1.2)
## scales 1.1.1 2020-05-11 [1] CRAN (R 4.1.0)
## sessioninfo 1.2.2 2021-12-06 [1] CRAN (R 4.1.0)
## stringi 1.7.6 2021-11-29 [1] CRAN (R 4.1.0)
## stringr * 1.4.0 2019-02-10 [1] CRAN (R 4.1.0)
## subplex 1.8 2022-04-12 [1] CRAN (R 4.1.2)
## testthat 3.1.3 2022-03-29 [1] CRAN (R 4.1.2)
## tibble * 3.1.6 2021-11-07 [1] CRAN (R 4.1.0)
## tidyr * 1.2.0 2022-02-01 [1] CRAN (R 4.1.2)
## tidyselect 1.1.2 2022-02-21 [1] CRAN (R 4.1.2)
## tidyverse * 1.3.1 2021-04-15 [1] CRAN (R 4.1.0)
## tzdb 0.3.0 2022-03-28 [1] CRAN (R 4.1.2)
## usethis 2.1.5 2021-12-09 [1] CRAN (R 4.1.0)
## utf8 1.2.2 2021-07-24 [1] CRAN (R 4.1.0)
## vctrs 0.4.0 2022-03-30 [1] CRAN (R 4.1.2)
## withr 2.5.0 2022-03-03 [1] CRAN (R 4.1.2)
## xfun 0.30 2022-03-02 [1] CRAN (R 4.1.2)
## xml2 1.3.3 2021-11-30 [1] CRAN (R 4.1.0)
## yaml 2.3.5 2022-02-21 [1] CRAN (R 4.1.2)
##
## [1] /Users/pals2/Library/R/x86_64/4.1/library
## [2] /Library/Frameworks/R.framework/Versions/4.1/Resources/library
##
## ──────────────────────────────────────────────────────────────────────────────