-
Notifications
You must be signed in to change notification settings - Fork 3
Expand file tree
/
Copy pathposterior_uncertainties.py
More file actions
133 lines (92 loc) · 3 KB
/
posterior_uncertainties.py
File metadata and controls
133 lines (92 loc) · 3 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
# ---
# jupyter:
# jupytext:
# text_representation:
# extension: .py
# format_name: light
# format_version: '1.5'
# jupytext_version: 1.15.2
# kernelspec:
# display_name: Python 3 (ipykernel)
# language: python
# name: python3
# ---
# # Posterior uncertainties
# In this tutorial, we will investigate the structure of the uncertainties
# returned by the probabilistic solvers. There is also the chance to compare
# filters and smoothers.
# +
"""Display the marginal uncertainties of filters and smoothers."""
# Set up the ODE
import jax
import jax.numpy as jnp
import matplotlib.pyplot as plt
from probdiffeq import ivpsolve, probdiffeq, taylor
# Fail this notebook on NaN detection (to catch those in the CI)
jax.config.update("jax_debug_nans", True)
@jax.jit
def vf(y, /, *, t): # noqa: ARG001
"""Evaluate the Lotka-Volterra vector field."""
y0, y1 = y[0], y[1]
y0_new = 0.5 * y0 - 0.05 * y0 * y1
y1_new = -0.5 * y1 + 0.05 * y0 * y1
return jnp.asarray([y0_new, y1_new])
t0 = 0.0
t1 = 2.0
u0 = jnp.asarray([20.0, 20.0])
# -
# Set up a solver.
#
# To all users: Try replacing the fixedpoint-smoother with a filter!
# +
tcoeffs = taylor.odejet_padded_scan(lambda y: vf(y, t=t0), (u0,), num=3)
init, ibm, ssm = probdiffeq.prior_wiener_integrated(tcoeffs, ssm_fact="blockdiag")
ts = probdiffeq.constraint_ode_ts1(vf, ssm=ssm)
strategy = probdiffeq.strategy_smoother_fixedpoint(ssm=ssm)
solver = probdiffeq.solver_mle(strategy=strategy, prior=ibm, constraint=ts, ssm=ssm)
errorest = probdiffeq.errorest_local_residual_cached(prior=ibm, ssm=ssm)
solve = ivpsolve.solve_adaptive_save_at(solver=solver, errorest=errorest)
# -
# Solve the ODE.
# +
ts = jnp.linspace(t0, t1, endpoint=True, num=50)
sol = jax.jit(solve)(init, save_at=ts, dt0=0.1, atol=1e-1, rtol=1e-1)
# -
# Plot the solution.
# +
fig, axes = plt.subplots(
nrows=3,
ncols=len(tcoeffs),
sharex="col",
tight_layout=True,
figsize=(len(sol.u.mean) * 2, 5),
)
for i, (u_i, std_i, ax_i) in enumerate(zip(sol.u.mean, sol.u.std, axes.T)):
# Set up titles and axis descriptions
if i == 0:
ax_i[0].set_title("State")
ax_i[0].set_ylabel("Prey")
ax_i[1].set_ylabel("Predators")
ax_i[2].set_ylabel("Std.-dev.")
elif i == 1:
ax_i[0].set_title(f"{i}st deriv.")
elif i == 2:
ax_i[0].set_title(f"{i}nd deriv.")
elif i == 3:
ax_i[0].set_title(f"{i}rd deriv.")
else:
ax_i[0].set_title(f"{i}th deriv.")
ax_i[-1].set_xlabel("Time")
for m, std, ax in zip(u_i.T, std_i.T, ax_i):
# Plot the mean
ax.plot(sol.t, m)
# Plot the standard deviation
lower, upper = m - 1.96 * std, m + 1.96 * std
ax.fill_between(sol.t, lower, upper, alpha=0.3)
ax.set_xlim((jnp.amin(ts), jnp.amax(ts)))
ax_i[2].semilogy(sol.t, std_i[:, 0], label="Prey")
ax_i[2].semilogy(sol.t, std_i[:, 1], label="Predators")
ax_i[2].legend(fontsize="x-small")
fig.align_ylabels()
plt.show()
# -