Implementing two changes:
This version fixes an error in the implementation of the non-central coefficients of the Fisher distributions as expressed in formula (9) and (10) in [1].
Using the notations from the article, theoretical formulas are (see [1], p.817, second column) :
$$ \left\{ \begin{aligned} \Delta_{10}^2(\mathbf{R}, \bar{\mathbf{N}}, \mathbf{C}^{\text{tip}}) & = \left\lVert(\mathbf{I} - \text{Proj}_{\mathbf{R}})\bar{\mathbf{N}}\right\lVert_{(\mathbf{C}^{\text{tip}})^{-1}}^2\\ \Delta_{20}^2(\mathbf{d}, \mathbf{R}, \bar{\mathbf{N}}, \mathbf{N}, \mathbf{C}^{\text{tip}}) & = \left\lVert(\mathbf{I} - \text{Proj}_{\mathbf{R}})\mathbf{N}\mathbf{d}\right\lVert_{(\mathbf{C}^{\text{tip}})^{-1}}^2 \;. \end{aligned} \right. $$where the projection is taken with respect to the metric defined by $(\mathbf{C}^{\text{tip}})^{-1}$, i.e. according to the formula:
$$ \text{Proj}_{\mathbf{M}} \mathbf{X} = \text{argmin}_{\mathbf{U}\in\text{Span}(\mathbf{M})} \left\lVert\mathbf{X} - \mathbf{U}\right\lVert_{(\mathbf{C}^{\text{tip}})^{-1}}^2 = \mathbf{M}(\mathbf{M}^T(\mathbf{C}^{\text{tip}})^{-1}\mathbf{M})^{-1}\mathbf{M}^T(\mathbf{C}^{\text{tip}})^{-1}\mathbf{X}\; . $$However, in the previous version of this notebook, the projection was taken with respect to the standard metric defined by the identity matrix.
Functions power10 and power20 of this version implement the correct formula.
Functions power10_OLD and power20_OLD are the previously used (erroneous) versions, kept for comparison.
The only impact this error has on the article is to change the estimated empirical power of the test $H_2$ versus $H_0$ for the sword index on the network with three hybridisation events from $0.47$ to $0.30$ (see below, and [1], p.812, first column, top).
This error has no impact on the theoretical power curves plotted in the article, that remain correct.
Julia 1.4.1.Initial release, along with the article:
[1] Paul Bastide, Claudia Solís-Lemus, Ricardo Kriebel, K William Sparks, Cécile Ané, Phylogenetic Comparative Methods on Phylogenetic Networks with Reticulations, Systematic Biology, Volume 67, Issue 5, September 2018, Pages 800–820 doi:10.1093/sysbio/syy033.
IJulia notebook compilation instructions¶This IJulia notebook (ipynb) file describes the
PCM analysis conducted on the fish genus Xiphophorus.
It was compiled into this html file using the Julia package
IJulia.
The interested reader can consult both:
.html file contains a summary of the most important commands
and results. It is self-sufficient..ipynb file can be opened with the IJulia notebook command, so that the commands included
in "chunks" of code can be executed. The default keystroke to execute a line of code from
the .ipynb file into the IJulia environment is Control + Enter.To launch the IJulia environement, the following commands can be executed.
(See the IJulia webpage for more information.)
# import Pkg
# Pkg.add("IJulia")
# using IJulia; notebook(detached=true)
PhyloNetworks package¶All the analyses presented here were done using functions from the
PhyloNetworks package
(Solís-Lemus et al, 2017).
A complete documentation is available
online.
To run the examples below, we need to install both PhyloNetworksand its
companion package PhyloPlots,
which contains all the network plotting features.
# import Pkg
# Pkg.add("PhyloNetworks")
# Pkg.add("PhyloPlots")
# Pkg.add("RCall")
after packages are installed, we can use them:
using PhyloNetworks
using PhyloPlots
using RCall
The specific versions of packages installed while runing this notebook are printed out at the end of this document, along with session informations.
In this section, we import and format the data from the online
supplementary files provided with this article. These commands all
assume that the associated files are in the same working directory as the
.ipynb file.
We fist read in the calibrated networks.
topologies = readMultiTopology("xiphophorus_networks_calibrated.tre");
tree = topologies[1];
net1 = topologies[2];
net3 = topologies[3];
All these networks share the same species. Next, we get the tip labels to use in matching and reducing the data.
taxa_net = tipLabels(tree);
We then import the morphology data.
using DataFrames
using CSV
dat = CSV.read("xiphophorus_morphology_Cui_etal_2013.csv",
types=[Union{Missing, String}, Union{Missing, Float64}, Union{Missing, Float64}],
copycols=true)
and delete the rows of data for the taxa not in the tree/networks.
for i in reverse(1:nrow(dat))
j = findfirst(taxa_net .== dat[!, :tipNames][i])
if isa(j, Nothing) # taxon not found in network
println("taxon ", dat[!, :tipNames][i], " (row $i) not found in network")
deleterows!(dat, i)
end
end
The first trait in the data we can study is the sword index. We conduct several PCM analyses using the three calibrated networks imported in the previous section.
We perform an ancestral state reconstruction assuming that the trait evolved like a BM on the network at hand. We start by taking a subset of the data with only the sword index.
dat_si = dat[:, [:sword_index, :tipNames]]
For each network, we plot both the predictions, and the 95% prediction intervals.
ASR_si_tree = ancestralStateReconstruction(dat_si, tree)
R"par(mar=c(.5,.5,.5,.5))" # to reduce margins of future plots
plot(tree, :RCall, nodeLabel=expectationsPlot(ASR_si_tree), useEdgeLength=true, xlim=[0,25]);
plot(tree, :RCall, nodeLabel=predintPlot(ASR_si_tree), useEdgeLength=true, xlim=[-3,26]);
ASR_si_net1 = ancestralStateReconstruction(dat_si, net1)
R"par(mar=c(.5,.5,.5,.5))"
plot(net1, :RCall, nodeLabel=expectationsPlot(ASR_si_net1), useEdgeLength=true, xlim=[0,25]);
plot(net1, :RCall, nodeLabel=predintPlot(ASR_si_net1), useEdgeLength=true, xlim=[-3,26]);
ASR_si_net3 = ancestralStateReconstruction(dat_si, net3)
R"par(mar=c(.5,.5,.5,.5))"
plot(net3, :RCall, nodeLabel=expectationsPlot(ASR_si_net3), useEdgeLength=true, xlim=[0,25]);
plot(net3, :RCall, nodeLabel=predintPlot(ASR_si_net3), useEdgeLength=true, xlim=[-3,26]);
using GLM
lambda_si_tree = phyloNetworklm(@formula(sword_index ~ 1), dat, tree, model="lambda")
round(lambda_estim(lambda_si_tree), digits = 2)
lambda_si_net1 = phyloNetworklm(@formula(sword_index ~ 1), dat, net1, model="lambda")
round(lambda_estim(lambda_si_net1), digits = 2)
lambda_si_net3 = phyloNetworklm(@formula(sword_index ~ 1), dat, net3, model="lambda")
round(lambda_estim(lambda_si_net3), digits = 2)
Transgressive evolution can only be tested for networks (h>0). The method relies on the construction of a regressor matrix that allows for a shift below each hybrid.
With a single reticulation event, we cannot test for heterogeneous effects of transgressive evolution. We can fit 2 models only: one with no transgressive evolution (H0), and another with a transgressive evolution effect after the 1 reticulation event (H1). Comparing these 2 models will be done with the F test for a (common) transgressive evolution effect.
We first construct the regressor matrix associated with hybrids.
regNet1 = regressorHybrid(net1)
This matrix has a column corresponding to each edge going out of an hybrid. It
indicates the tips descending from each, taking the possible partial
descendances introduced by the network into account. It also has an extra
column labelled sum, that we will discuss afterward. Here, only the edge
with number 34 is below a reticulation, as we can see from plotting the network
using edge numbers (where we do not use the branch lengths, so that the
topology appears more clearly). The associated predictor is labelled shift_34
(see above).
R"par(mar=c(.5,.5,.5,.5))"
plot(net1, :RCall, showEdgeNumber=true, xlim=[0,15]);
We can then join the regressor matrix with the data, in order to
use the phyloNetworklm function. We call this function twice.
First, we fit the data against a simple intercept (H0).
Then, we fit the data against the regressor induced by the shift after the
hybrid (H1).
Once this is done, we simply call the GLM function ftest to perform
the F test, including the models from the most to the least complex.
dat1 = join(dat, regNet1, on = :tipNames);
fit_net1_0 = phyloNetworklm(@formula(sword_index ~ 1), dat1, net1)
fit_net1_1 = phyloNetworklm(@formula(sword_index ~ shift_34), dat1, net1)
ftest(fit_net1_1, fit_net1_0)
Note that, for conventional reasons, the ftest function takes the
most complex model as the first one. This means that, in the table of
results, the models are actually named in a reverse order, so that
"Model 2" is actually our model under H0, and "Model 1" the one under
H1.
Here with 3 reticulations, we can consider 3 models: no transgressive evolution (H0), common transgressive evolution (H1), and heterogeneous transgressive evolution (H2). We can run F tests to compare H1 vs H0, and H2 vs H1.
We construct the regressor matrix the same way as before.
regNet3 = regressorHybrid(net3)
Here, 3 columns are created, shift_24, shift_37 and shift_45, associated
with the 3 hybridization events occurring in the network. The extra column
sum is just the sum of these 3 predictors. It is useful to perform the fit of
model H1, where we only consider a homogeneous transgressive evolution effect.
For the third fit (H2), we need to consider all three predictors independently.
The results of the F tests can hence be obtained as follow.
dat3 = join(dat, regNet3, on = :tipNames);
fit_net3_0 = phyloNetworklm(@formula(sword_index ~ 1), dat3, net3)
fit_net3_1 = phyloNetworklm(@formula(sword_index ~ sum), dat3, net3)
fit_net3_2 = phyloNetworklm(@formula(sword_index ~ shift_24 + shift_37 + shift_45),
dat3, net3)
ftest(fit_net3_2, fit_net3_1, fit_net3_0)
As above, beware that the numbering of the models in the table of result (most to least complex) is reversed compared to our numbering (least to most complex). Hence, here, "Model 3" is actually our model under H0, "Model 2" the one under H1, and "Model 1" the one under H2.
We can comute the power associated with the tests above assuming that the network, variance parameter and size of the effects are correctly estimated.
We first define two functions computing the power of the homogeneous vs no
effect (power10) and heterogeneous vs no effect (power20),
when the p-value is called significant if lower than α=0.05
(but change the "level" argument below to change α).
using Distributions
using LinearAlgebra
function power10_OLD(b::Real, sigma2::Real, net::HybridNetwork; level=0.95)
# Characteristics of the network
n = length(tipLabels(net))
V = sharedPathMatrix(net)[:Tips]
dfr_hybrid = regressorHybrid(net)
# Compute the non-central coefficient
X = ones(n)
Z = convert(Array, dfr_hybrid[!, :sum])
ncc10 = (I - 1/n .* X*X') * Z ## WRONG PROJECTION HERE
ncc10 = ncc10' * V^(-1) * ncc10
ncc10 = b^2/sigma2 * ncc10[1]
# Compute power
d1 = 2 - 1
d2 = n - 2
q95 = quantile(FDist(d1, d2), level)
pow95 = ccdf(NoncentralF(d1, d2, ncc10), q95)
end
function power10(b::Real, sigma2::Real, net::HybridNetwork; level=0.95)
# Characteristics of the network
n = length(tipLabels(net))
V = sharedPathMatrix(net)[:Tips]
dfr_hybrid = regressorHybrid(net)
# Compute the non-central coefficient
X = ones(n)
Z = convert(Array, dfr_hybrid[!, :sum])
Vinv = V^(-1)
ncc10 = (I - X*(X'*Vinv*X)^(-1)*X'*Vinv) * Z ## CORRECT PROJECTION HERE
ncc10 = ncc10' * Vinv * ncc10
ncc10 = b^2/sigma2 * ncc10[1]
# Compute power
d1 = 2 - 1
d2 = n - 2
q95 = quantile(FDist(d1, d2), level)
pow95 = ccdf(NoncentralF(d1, d2, ncc10), q95)
end
function power20_OLD(B::Vector, sigma2::Real, net::HybridNetwork; level=0.95)
# Characteristics of the network
n = length(tipLabels(net))
V = sharedPathMatrix(net)[:Tips]
dfr_hybrid = regressorHybrid(net)
# Regression matrices
X = [ones(n) convert(Array, dfr_hybrid[!, :sum])]
Z = convert(Array, dfr_hybrid[:, filter(x -> !(x in [:sum, :tipNames]), names(dfr_hybrid))])
nh = size(Z,2)
# noncentral coef
ncc20 = (I - X*(X'*X)^(-1)*X') * Z * B ## WRONG PROJECTION HERE
ncc20 = ncc20' * V^(-1) * ncc20
ncc20 = 1/sigma2 * ncc20[1]
# power
d1 = nh - 1
d2 = n - nh - 1
q95 = quantile(FDist(d1, d2), level)
pow95 = ccdf(NoncentralF(d1, d2, ncc20), q95)
end
function power20(B::Vector, sigma2::Real, net::HybridNetwork; level=0.95)
# Characteristics of the network
n = length(tipLabels(net))
V = sharedPathMatrix(net)[:Tips]
dfr_hybrid = regressorHybrid(net)
# Regression matrices
X = [ones(n) convert(Array, dfr_hybrid[!, :sum])]
Z = convert(Array, dfr_hybrid[:, filter(x -> !(x in [:sum, :tipNames]), names(dfr_hybrid))])
nh = size(Z,2)
# noncentral coef
Vinv = V^(-1)
ncc20 = (I - X*(X'*Vinv*X)^(-1)*X'*Vinv) * Z * B ## CORRECT PROJECTION HERE
ncc20 = ncc20' * Vinv * ncc20
ncc20 = 1/sigma2 * ncc20[1]
# power
d1 = nh - 1
d2 = n - nh - 1
q95 = quantile(FDist(d1, d2), level)
pow95 = ccdf(NoncentralF(d1, d2, ncc20), q95)
end
power10(coef(fit_net1_1)[2], sigma2_estim(fit_net1_1), net1)
power10_OLD(coef(fit_net1_1)[2], sigma2_estim(fit_net1_1), net1)
power10(coef(fit_net3_1)[2], sigma2_estim(fit_net3_1), net3)
power10_OLD(coef(fit_net3_1)[2], sigma2_estim(fit_net3_1), net3)
power20(coef(fit_net3_2)[2:4], sigma2_estim(fit_net3_2), net3)
power20_OLD(coef(fit_net3_2)[2:4], sigma2_estim(fit_net3_2), net3)
The power to detect homogeneous effect is very small (around 10%),
when using parameters estimated from the sword index data.
This is consistent with the graph presented in the main text,
as one need to use the normalized coefficient
$$\sqrt{\frac{b^2}{\sigma^2 h}}$$
where h is the total heigh of the tree.
Here, for net1, this coefficient is about 0.2 (computed below).
sqrt(coef(fit_net3_1)[2]^2/sigma2_estim(fit_net3_1)/maximum(getNodeAges(net3)))
However, the power to detect heterogeneous effects is a bit better, around 47%.
We reproduce these analyses for the second trait, the female preference.
Some taxa present in the network are missing data for female preference. In all this section, we do not remove the taxa with missing data, from the network or from the morphology data frame.
That enables us to:
As above, we assume that the trait evolved like a BM on the network at hand. We start by taking a sub-set of the data with only the female preference.
dat_fp = dat[:, [:preference, :tipNames]]
For each network, we plot both the predictions, and the 95% prediction intervals.
ASR_fp_tree = ancestralStateReconstruction(dat_fp, tree)
R"par(mar=c(.5,.5,.5,.5))"
plot(tree, :RCall, nodeLabel=expectationsPlot(ASR_fp_tree), useEdgeLength=true, xlim=[0,25]);
plot(tree, :RCall, nodeLabel=predintPlot(ASR_fp_tree), useEdgeLength=true, xlim=[-3,26]);
Note that imputed trait values are also plotted with their prediction intervals.
ASR_fp_net1 = ancestralStateReconstruction(dat_fp, net1)
R"par(mar=c(.5,.5,.5,.5))"
plot(net1, :RCall, nodeLabel=expectationsPlot(ASR_fp_net1), useEdgeLength=true, xlim=[0,25]);
plot(net1, :RCall, nodeLabel=predintPlot(ASR_fp_net1), useEdgeLength=true, xlim=[-3,26]);
ASR_fp_net3 = ancestralStateReconstruction(dat_fp, net3)
R"par(mar=c(.5,.5,.5,.5))"
plot(net3, :RCall, nodeLabel=expectationsPlot(ASR_fp_net3), useEdgeLength=true, xlim=[0,25]);
plot(net3, :RCall, nodeLabel=predintPlot(ASR_fp_net3), useEdgeLength=true, xlim=[-3,26]);
lambda_fp_tree = phyloNetworklm(@formula(preference ~ 1), dat, tree, model="lambda")
round(lambda_estim(lambda_fp_tree), digits = 2)
lambda_fp_net1 = phyloNetworklm(@formula(preference ~ 1), dat, net1, model="lambda")
round(lambda_estim(lambda_fp_net1), digits = 2)
lambda_fp_net3 = phyloNetworklm(@formula(preference ~ 1), dat, net3, model="lambda")
round(lambda_estim(lambda_fp_net3), digits = 2)
We apply the same methods to study transgressive evolution on both networks. The approach relies on the construction of the regressor matrix, and then on a fit of the different models, on which a F test can be performed.
Here, the regression matrices are the same as above, as the networks are the same. We can hence use the ones we already computed. We just recall the way to do it here for clarity, but will not show the regressor matrices.
regNet1 = regressorHybrid(net1)
dat1 = join(dat, regNet1, on = :tipNames);
fit_net1_0 = phyloNetworklm(@formula(preference ~ 1), dat1, net1)
fit_net1_1 = phyloNetworklm(@formula(preference ~ shift_34), dat1, net1)
ftest(fit_net1_1, fit_net1_0)
regNet3 = regressorHybrid(net3)
dat3 = join(dat, regNet3, on = :tipNames);
fit_net3_0 = phyloNetworklm(@formula(preference ~ 1), dat3, net3)
fit_net3_1 = phyloNetworklm(@formula(preference ~ sum), dat3, net3)
fit_net3_2 = phyloNetworklm(@formula(preference ~ shift_24 + shift_37 + shift_45),
dat3, net3)
ftest(fit_net3_2, fit_net3_1, fit_net3_0)
Here, we notice that the F test indicates a non-significant homogeneous
effect (H1 vs H0, i.e. Model 2 vs Model 3 in the table), but a significant heterogeneous effect
(H2 vs H1, i.e. Model 1 vs Model 2 in the table).
It might be informative to test directly the full model against the null (H2 vs H0).
This can be done easily applying the ftest function to the corresponding fits.
ftest(fit_net3_2, fit_net3_0)
As above, we can compute the power associated with the tests above assuming that the variance paremter and size of the effects are correctly estimated.
power10(coef(fit_net1_1)[2], sigma2_estim(fit_net1_1), net1)
power10_OLD(coef(fit_net1_1)[2], sigma2_estim(fit_net1_1), net1)
power10(coef(fit_net3_1)[2], sigma2_estim(fit_net3_1), net3)
power10_OLD(coef(fit_net3_1)[2], sigma2_estim(fit_net3_1), net3)
power20(coef(fit_net3_2)[2:4], sigma2_estim(fit_net3_2), net3)
power20_OLD(coef(fit_net3_2)[2:4], sigma2_estim(fit_net3_2), net3)
Here, the power to detect effect is much higher, going up to almost 1 for heterogeneous effects. This is consistent with the effect of transgressive evolution being higher.
sqrt(coef(fit_net3_1)[2]^2/sigma2_estim(fit_net3_1)/maximum(getNodeAges(net3)))
We reproduce here the same analysis, but this time removing missing data from the network. Pruning taxa from networks is done below.
"""
reduceNet!(phynet)
function to match data to a network, to
- delete tips in network if they lack "preference" data, and
- delete rows in data "dat" if not in network.
notes:
- reticulations reducing to a simple loop are simplified (removed)
- object "dat" is used (and modified) as an external variable,
it needs to be defined beforehand.
"""
function reduceNet!(phynet)
taxa = tipLabels(phynet)
for i in 1:phynet.numTaxa
j = findfirst(dat[:tipNames], taxa[i])
if j==0 # taxon not found in trait data set
#println("tipNames ", taxa[i], " not found in trait data")
deleteleaf!(phynet, taxa[i])
continue # skips next part: where j needs to be >0
end
if ismissing(dat[:preference][j]) # taxon found but has missing trait value
#println("tipNames ", taxa[i], " lacks trait value")
deleteleaf!(phynet, taxa[i])
deleterows!(dat, j)
end
end
end
We now apply our function to each network:
reduceNet!(tree)
reduceNet!(net1)
reduceNet!(net3)
Once we modified the data and networks to remove taxa with missing preference data, the analysis are almost exactly the same as in the previous section, only the results differ.
As above, we assume that the trait evolved like a BM on the network at hand. We start by taking a subset of the data with only the female preference.
dat_fp = dat[:, [:preference, :tipNames]]
For each network, we plot both the predictions, and the 95% prediction intervals.
ASR_fp_tree = ancestralStateReconstruction(dat_fp, tree)
R"par(mar=c(.5,.5,.5,.5))"
plot(tree, :RCall, nodeLabel=expectationsPlot(ASR_fp_tree), useEdgeLength=true, xlim=[0,25]);
plot(tree, :RCall, nodeLabel=predintPlot(ASR_fp_tree), useEdgeLength=true, xlim=[-3,26]);
ASR_fp_net1 = ancestralStateReconstruction(dat_fp, net1)
R"par(mar=c(.5,.5,.5,.5))"
plot(net1, :RCall, nodeLabel=expectationsPlot(ASR_fp_net1), useEdgeLength=true, xlim=[0,25]);
plot(net1, :RCall, nodeLabel=predintPlot(ASR_fp_net1), useEdgeLength=true, xlim=[-3,26]);
ASR_fp_net3 = ancestralStateReconstruction(dat_fp, net3)
R"par(mar=c(.5,.5,.5,.5))"
plot(net3, :RCall, nodeLabel=expectationsPlot(ASR_fp_net3), useEdgeLength=true, xlim=[0,25]);
plot(net3, :RCall, nodeLabel=predintPlot(ASR_fp_net3), useEdgeLength=true, xlim=[-3,26]);
lambda_fp_pruned_tree = phyloNetworklm(@formula(preference ~ 1), dat, tree, model="lambda")
round(lambda_estim(lambda_fp_pruned_tree), digits = 2)
lambda_fp_pruned_net1 = phyloNetworklm(@formula(preference ~ 1), dat, net1, model="lambda")
round(lambda_estim(lambda_fp_pruned_net1), digits = 2)
lambda_fp_pruned_net3 = phyloNetworklm(@formula(preference ~ 1), dat, net3, model="lambda")
round(lambda_estim(lambda_fp_pruned_net3), digits = 2)
regNet1 = regressorHybrid(net1)
dat1 = join(dat, regNet1, on = :tipNames);
fit_net1_0 = phyloNetworklm(@formula(preference ~ 1), dat1, net1);
fit_net1_1 = phyloNetworklm(@formula(preference ~ shift_34), dat1, net1);
ftest(fit_net1_1, fit_net1_0)
When pruning the second network, one of the hybridization
event was reduced to a simple loop over itself, and hence was dropped
from the network. As a result, net3 has now only two hybridization events,
as we can see on the following plot.
R"par(mar=c(.5,.5,.5,.5))"
plot(net3, :RCall, showEdgeNumber=true, xlim=[0,15]);
From this topology, we can conclude that only two predictors will be included in the regression matrix, associated with edges 24 and 37. The analysis can hence be conducted the following way.
regNet3 = regressorHybrid(net3)
dat3 = join(dat, regNet3, on = :tipNames);
fit_net3_0 = phyloNetworklm(@formula(preference ~ 1), dat3, net3)
fit_net3_1 = phyloNetworklm(@formula(preference ~ sum), dat3, net3)
fit_net3_2 = phyloNetworklm(@formula(preference ~ shift_24 + shift_37), dat3, net3)
ftest(fit_net3_2, fit_net3_1, fit_net3_0)
Here, contrary to the results above, we get a significant homogeneous effect, but no evidence for heterogeneous effects.
This hypothesis was originally tested by Cui et al. 2013. We can plot the data first:
using Gadfly
Gadfly.plot(dat, x="sword_index", y="preference", Guide.ylabel("female preference"))
phyloNetworklm(@formula(sword_index ~ preference), dat, tree)
phyloNetworklm(@formula(sword_index ~ preference), dat, net1)
phyloNetworklm(@formula(sword_index ~ preference), dat, net3)
On all networks, we see a positive relationship with slopes of 0.5878, 0.5529 and 0.5563. But this relationship is not significant (p=0.096, 0.110, 0.106).
import Dates
Dates.now()
versioninfo()
import Pkg
Pkg.status()