Figure 1 describes the procedure we will pursue:
define the graph and add blocks (sections A, B and C in the Figure)
simulate the model, instead we could also use actual data (section D in Figure)
compute the cross spectral density
setup the DCM
estimate parameters
plot the results
Learning goals
perform the entire workflow of an spDCM analysis.
use observer Blox to simulate experimental measurements.
using Neuroblox
using LinearAlgebra
using StochasticDiffEq
using DataFrames
using OrderedCollections
using CairoMakie
using ModelingToolkit
using Random
We will define a model of 3 regions. This means first of all to define a graph. To this graph we will add three linear neuronal mass models which constitute the (hidden) neuronal dynamics. These constitute three nodes of the graph. Next we will also need some input that stimulates the activity, we use simple Ornstein-Uhlenbeck blocks to create stochastic inputs. One per region. We want to simulate fMRI signals thus we will need to also add a BalloonModel per region. Note that the Ornstein-Uhlenbeck block will feed into the linear neural mass which in turn will feed into the BalloonModel blox. This needs to be represented by the way we define the edges.
Random.seed!(17) # set seed for reproducibility
nr = 3 # number of regions
g = MetaDiGraph()
regions = []; # list of neural mass blocks to then connect them to each other with an adjacency matrix `A_true`
Now add the different blocks to each region and connect the blocks within each region. For convenience we use a for loop since the type of blocks belonging to a each region repeat over regions but you could also approach building the system the same way as was shown in previous tutorials:
for i = 1:nr
region = LinearNeuralMass(;name=Symbol("r$(i)₊lm"))
push!(regions, region) # store neural mass model in list. We need this list below. If you haven't seen the Julia command `push!` before [see here](http://jlhub.com/julia/manual/en/function/push-exclamation).
input = OUBlox(;name=Symbol("r$(i)₊ou"), σ=0.1) ## add Ornstein-Uhlenbeck as noisy input to the current region
add_edge!(g, input => region, weight=1/16)
measurement = BalloonModel(;name=Symbol("r$(i)₊bm")) ## simulate fMRI signal with BalloonModel which includes the BOLD signal on top of the balloon model dynamics
add_edge!(g, region => measurement, weight=1.0)
end
Note that weight=1/16
in the connection between the OU process and the Neural Mass Blox is taken from SPM12. This stabilizes the balloon model simulation. Alternatively the noise of the Ornstein-Uhlenbeck block or the weight of the edge connecting neuronal activity and balloon model could be reduced to guarantee numerical stability. Next we define the between-region connectivity matrix and connect regions; we use the same matrix as is defined in [3]
A_true = [[-0.5 -2 0]; [0.4 -0.5 -0.3]; [0 0.2 -0.5]]
for idx in CartesianIndices(A_true)
add_edge!(g, regions[idx[1]] => regions[idx[2]], weight=A_true[idx[1], idx[2]])
end
finally we compose the simulation model
@named simmodel = system_from_graph(g);
setup simulation of the model, time in seconds
tspan = (0.0, 512.0)
prob = SDEProblem(simmodel, [], tspan)
dt = 2 # 2 seconds (units are seconds) as measurement interval for fMRI
sol = solve(prob, ImplicitRKMil(), saveat=dt);
we now want to extract all the variables in our model which carry the tag "measurement". For this purpose we can use the Neuroblox function get_idx_tagged_vars
the observable quantity in our model is the BOLD signal, the variable of the Blox BalloonModel
that represents the BOLD signal is tagged with "measurement" tag. other tags that are defined are "input" which denotes variables representing a stimulus, like for instance an OUBlox
.
idx_m = get_idx_tagged_vars(simmodel, "measurement") # get index of bold signal
3-element Vector{Int64}:
19
20
21
plot bold signal time series
f = Figure()
ax = Axis(f[1, 1],
title = "fMRI time series",
xlabel = "Time [ms]",
ylabel = "BOLD",
)
lines!(ax, sol, idxs=idx_m)
f
We note that the initial spike is not meaningful and a result of the equilibration of the stochastic process thus we remove it.
dfsol = DataFrame(sol[ceil(Int, 101/dt):end]);
data = Matrix(dfsol[:, idx_m]);
We compute the cross-spectral density by fitting a linear model of order p
and then compute the csd analytically from the parameters of the multivariate autoregressive model
p = 8
mar = mar_ml(data, p) # maximum likelihood estimation of the MAR coefficients and noise covariance matrix
ns = size(data, 1)
freq = range(min(128, ns*dt)^-1, max(8, 2*dt)^-1, 32)
csd = mar2csd(mar, freq, dt^-1);
Now plot the cross-spectrum:
fig = Figure(size=(1200, 800))
grid = fig[1, 1] = GridLayout()
for i = 1:nr
for j = 1:nr
ax = Axis(grid[i, j])
lines!(ax, freq, real.(csd[:, i, j]))
end
end
fig
We will now assemble a new model that is used for fitting the previous simulations. This procedure is similar to before with the difference that we will define global parameters and use tags such as [tunable=false/true] to define which parameters we will want to estimate. Note that parameters are tunable by default.
g = MetaDiGraph()
regions = []; # list of neural mass blocks to then connect them to each other with an adjacency matrix `A`
Note that parameters are typically defined within a Blox and thus not immediately visible to the user. Since we want some parameters to be shared across several regions we define them outside of the regions. For this purpose use the ModelingToolkit macro @parameters
which is used to define symbolic parameters for models. Note that we can set the tunable flag right away thereby defining whether we will include this parameter in the optimization procedure or rather keep it fixed to its predefined value.
@parameters lnκ=0.0 [tunable=false] lnϵ=0.0 [tunable=false] lnτ=0.0 [tunable=false]; # lnκ: decay parameter for hemodynamics; lnϵ: ratio of intra- to extra-vascular components, lnτ: transit time scale
@parameters C=1/16 [tunable=false]; # note that C=1/16 is taken from SPM12 and stabilizes the balloon model simulation. See also comment above.
We now define a similar model as above for the simulation but instead of using an actual stimulus Blox we here add ExternalInput which represents a simple linear external input that is not specified any further. We simply say that our model gets some input with a proportional factor . This is mostly only to make sure that our results are consistent with those produced by SPM
for i = 1:nr
region = LinearNeuralMass(;name=Symbol("r$(i)₊lm"))
push!(regions, region)
input = ExternalInput(;name=Symbol("r$(i)₊ei"))
add_edge!(g, input => region, weight=C)
measurement = BalloonModel(;name=Symbol("r$(i)₊bm"), lnτ=lnτ, lnκ=lnκ, lnϵ=lnϵ) ## assume fMRI signal and model them with a BalloonModel
add_edge!(g, region => measurement, weight=1.0)
end
Here we define the prior expectation values of the effective connectivity matrix we wish to infer:
A_prior = 0.01*randn(nr, nr)
A_prior -= diagm(diag(A_prior)) # remove the diagonal
3×3 Matrix{Float64}:
0.0 0.00345049 -0.0122152
0.0148102 0.0 0.00437472
-0.00716553 0.00429191 0.0
Since we want to optimize these weights we turn them into symbolic parameters: Add the symbolic weights to the edges and connect regions.
A = []
for (i, a) in enumerate(vec(A_prior))
symb = Symbol("A$(i)")
push!(A, only(@parameters $symb = a))
end
for (i, idx) in enumerate(CartesianIndices(A_prior))
if idx[1] == idx[2]
add_edge!(g, regions[idx[1]] => regions[idx[2]], weight=-exp(A[i])/2) ## -exp(A[i])/2: treatement of diagonal elements in SPM12 to make diagonal dominance (see Gershgorin Theorem) more likely but it is not guaranteed
else
add_edge!(g, regions[idx[2]] => regions[idx[1]], weight=A[i])
end
end
# Avoid simplification of the model in order to be able to exclude some parameters from fitting
@named fitmodel = system_from_graph(g, simplify=false);
With the Neuroblox function changetune
we can provide a dictionary of parameters whose tunable flag should be changed, for instance set to false to exclude them from the optimization procedure. Assume we want to exclude the connections that were set to zero in the simulation:
untune = Dict(A[3] => false, A[7] => false)
fitmodel = changetune(fitmodel, untune) # A[3] and A[7] were set to 0 in the simulation
fitmodel = structural_simplify(fitmodel, split=false) # and now simplify the euqations; the `split` parameter is necessary for some ModelingToolkit peculiarities and will soon be removed. So don't lose time with it ;)
Model fitmodel:
Equations (21):
21 standard: see equations(fitmodel)
Unknowns (21): see unknowns(fitmodel)
r1₊lm₊x(t) [defaults to 0.0]
r1₊bm₊s(t) [defaults to 1.0]
r1₊bm₊lnu(t) [defaults to 0.0]
r1₊bm₊lnν(t) [defaults to 0.0]
r1₊bm₊lnq(t) [defaults to 0.0]
r2₊lm₊x(t) [defaults to 0.0]
r2₊bm₊s(t) [defaults to 1.0]
r2₊bm₊lnu(t) [defaults to 0.0]
r2₊bm₊lnν(t) [defaults to 0.0]
r2₊bm₊lnq(t) [defaults to 0.0]
r3₊lm₊x(t) [defaults to 0.0]
r3₊bm₊s(t) [defaults to 1.0]
r3₊bm₊lnu(t) [defaults to 0.0]
r3₊bm₊lnν(t) [defaults to 0.0]
r3₊bm₊lnq(t) [defaults to 0.0]
r1₊ei₊u(t): ext_input
r2₊ei₊u(t): ext_input
r3₊ei₊u(t): ext_input
r1₊bm₊bold(t) [defaults to 0.0]: measurement
r2₊bm₊bold(t) [defaults to 0.0]: measurement
r3₊bm₊bold(t) [defaults to 0.0]: measurement
Parameters (16): see parameters(fitmodel)
C [defaults to 0.0625]
A1 [defaults to 0.0]
w_r1₊lm_r1₊bm [defaults to 1.0]
A2 [defaults to 0.0148102]
A3 [defaults to -0.00716553]
A4 [defaults to 0.00345049]
A5 [defaults to 0.0]
w_r2₊lm_r2₊bm [defaults to 1.0]
A6 [defaults to 0.00429191]
A7 [defaults to -0.0122152]
A8 [defaults to 0.00437472]
A9 [defaults to 0.0]
w_r3₊lm_r3₊bm [defaults to 1.0]
lnκ [defaults to 0.0]
lnτ [defaults to 0.0]
lnϵ [defaults to 0.0]
Observed (6): see observed(fitmodel)
max_iter = 128; # maximum number of iterations
# attribute initial conditions or default values to dynamic states of our model
sts, _ = get_dynamic_states(fitmodel);
the following step is needed if the model's Jacobian would give degenerate eigenvalues when expanded around the fixed point 0 (which is the default expansion). We simply add small random values to avoid this degeneracy:
perturbedfp = Dict(sts .=> abs.(0.001*rand(length(sts)))) # slight noise to avoid issues with Automatic Differentiation.
Dict{SymbolicUtils.BasicSymbolic{Real}, Float64} with 15 entries:
r3₊lm₊x(t) => 0.000489568
r3₊bm₊lnq(t) => 0.000198145
r2₊bm₊lnq(t) => 0.000210419
r1₊bm₊lnν(t) => 0.000541726
r2₊bm₊s(t) => 0.000397204
r3₊bm₊lnu(t) => 0.000169086
r2₊bm₊lnν(t) => 0.000646358
r1₊bm₊s(t) => 0.000138322
r3₊bm₊lnν(t) => 0.000570881
r1₊lm₊x(t) => 4.49163e-6
r1₊bm₊lnu(t) => 0.000896385
r3₊bm₊s(t) => 0.000388604
r2₊lm₊x(t) => 0.000928402
r2₊bm₊lnu(t) => 0.000747008
r1₊bm₊lnq(t) => 0.000808084
For convenience we can use the default prior function to use standardized prior values as given in SPM:
pmean, pcovariance, indices = defaultprior(fitmodel, nr)
priors = (μθ_pr = pmean,
Σθ_pr = pcovariance
);
Setup hyper parameter prior as well:
hyperpriors = Dict(:Πλ_pr => 128.0*ones(1, 1), # prior metaparameter precision, needs to be a matrix
:μλ_pr => [8.0] # prior metaparameter mean, needs to be a vector
);
To compute the cross spectral densities we need to provide the sampling interval of the time series, the frequency axis and the order of the multivariate autoregressive model:
csdsetup = (mar_order = p, freq = freq, dt = dt);
earlier we used the function get_idx_tagged_vars
to get the indices of tagged variables. Here we don't want to get the indices but rather the symbolic variable names themselves in order to get the correct columns of the dataframe of the simulation that correspond to the BOLD signal or measurement:
_, s_bold = get_eqidx_tagged_vars(fitmodel, "measurement"); # get bold signal variables
Prepare the DCM. This function will setup the computation of the Dynamic Causal Model. The last parameter specifies that we are using fMRI time series (as opposed to LFPs, which is the other modality that is currently available in Neuroblox).
(state, setup) = setup_sDCM(dfsol[:, String.(Symbol.(s_bold))], fitmodel, perturbedfp, csdsetup, priors, hyperpriors, indices, pmean, "fMRI");
We are now ready to run the optimization procedure! That is we loop over runsDCMiteration! which will alter state
after each optimization iteration. It essentially computes the Variational Laplace estimation of expectation and variance of the tunable parameters.
for iter in 1:max_iter
state.iter = iter
run_sDCM_iteration!(state, setup)
print("iteration: ", iter, " - F:", state.F[end], " - dF predicted:", state.dF[end], "\n")
if iter >= 4
criterion = state.dF[end-3:end] .< setup.tolerance
if all(criterion)
print("convergence\n")
break
end
end
end
iteration: 1 - F:-3568.7760918985473 - dF predicted:1151.645766454626
iteration: 2 - F:-2848.2213182820224 - dF predicted:570.6511079776099
iteration: 3 - F:-2353.1242235135683 - dF predicted:524.3646402335777
iteration: 4 - F:-1903.0539222090856 - dF predicted:478.7683995366398
iteration: 5 - F:-1485.7828414862329 - dF predicted:453.24897986787704
iteration: 6 - F:-1097.437817610763 - dF predicted:365.0822004373764
iteration: 7 - F:-804.3167412039838 - dF predicted:210.2013606453641
iteration: 8 - F:-644.0670887998349 - dF predicted:100.06125194194902
iteration: 9 - F:-564.1829538619753 - dF predicted:53.19735621850561
iteration: 10 - F:-522.2850500695243 - dF predicted:24.073777898263447
iteration: 11 - F:-504.33225663624694 - dF predicted:8.902065051716994
iteration: 12 - F:-497.4118392749408 - dF predicted:3.901165216297128
iteration: 13 - F:-494.19989921140956 - dF predicted:2.2364255974822576
iteration: 14 - F:-492.321101184751 - dF predicted:1.429503973813221
iteration: 15 - F:-491.1320145742268 - dF predicted:0.8876270589104192
iteration: 16 - F:-490.4043700125305 - dF predicted:0.5419600356351018
iteration: 17 - F:-489.95846375941073 - dF predicted:0.32462393588823557
iteration: 18 - F:-489.7081937123439 - dF predicted:0.1953923146434604
iteration: 19 - F:-489.5533654142472 - dF predicted:0.135656721526179
iteration: 20 - F:-489.50792749644734 - dF predicted:0.18761426800246586
iteration: 21 - F:-489.498524195646 - dF predicted:0.6009586126933194
iteration: 22 - F:-489.498524195646 - dF predicted:0.002476062697706699
iteration: 23 - F:-489.4966992190876 - dF predicted:0.0037885150812212915
iteration: 24 - F:-489.493464681141 - dF predicted:0.005932986591946685
iteration: 25 - F:-489.48824801595504 - dF predicted:0.009058419422460832
convergence
Note that the output F
is the free energy at each iteration step and dF
is the predicted change of free energy at each step which approximates the actual free energy change and is used as stopping criterion by requiring that it does not excede the tolerance
level for 4 consecutive times.
Free energy is the objective function of the optimization scheme of spectral DCM. Note that in the machine learning literature this it is called Evidence Lower Bound (ELBO). Plot the free energy evolution over optimization iterations to see how the algorithm converges towards a (potentially local) optimum:
f1 = freeenergy(state)
Plot the estimated posterior of the effective connectivity and compare that to the true parameter values. Bar height are the posterior mean and error bars are the standard deviation of the posterior.
f2 = ecbarplot(state, setup, A_true)
Explore susceptibility with respect to noise. Run the script again with a different random seed and observe how the results change. Given that we didn’t change any parameters of the ground truth, what is your take on parameter inference with this setup? How reliable is model selection based on free energy (compare the different free energies of the models and their respective parameter values with ground truth)?
Averaging over patients. Now repeat the simulation and inference again with 10 different seeds (you can just remove the number in Random.seed() to use current time stamps as seeds) and store the results of each run, such that you get several models based on the same ground truth but with different instances of the noise. Can you extract the average value of the effective connectivity from the ensemble?
Changing neuronal dynamics model. Now change the model and test the whole procedure with a different underlying neuronal mass model, for instance the Jansen-Rit model. Note that there are no default priors for the Jansen-Rit model, you will have to provide priors for the extra parameters or remove them from the optimization procedure by setting their tunable to false.