Skip to content

Commit 70218bf

Browse files
ENH: add exponentiation of a covariance function with a scalar (pymc-devs#3852)
* ENH: add exponentiation with a scalar * fix the scalar condition * add examples in notebook
1 parent c34ae3f commit 70218bf

File tree

4 files changed

+271
-37
lines changed

4 files changed

+271
-37
lines changed

RELEASE-NOTES.md

+1
Original file line numberDiff line numberDiff line change
@@ -30,6 +30,7 @@
3030
- Distinguish between `Data` and `Deterministic` variables when graphing models with graphviz. PR [#3491](https://github.com/pymc-devs/pymc3/pull/3491).
3131
- Sequential Monte Carlo - Approximate Bayesian Computation step method is now available. The implementation is in an experimental stage and will be further improved.
3232
- Added `Matern12` covariance function for Gaussian processes. This is the Matern kernel with nu=1/2.
33+
- Added exponentiation of a covariance function with a scalar. See PR[#3852](https://github.com/pymc-devs/pymc3/pull/3852)
3334
- Progressbar reports number of divergences in real time, when available [#3547](https://github.com/pymc-devs/pymc3/pull/3547).
3435
- Sampling from variational approximation now allows for alternative trace backends [#3550].
3536
- Infix `@` operator now works with random variables and deterministics [#3619](https://github.com/pymc-devs/pymc3/pull/3619).

docs/source/notebooks/GP-MeansAndCovs.ipynb

+184-37
Large diffs are not rendered by default.

pymc3/gp/cov.py

+31
Original file line numberDiff line numberDiff line change
@@ -13,9 +13,11 @@
1313
# limitations under the License.
1414

1515
import numpy as np
16+
import theano
1617
import theano.tensor as tt
1718
from functools import reduce
1819
from operator import mul, add
20+
from numbers import Number
1921

2022
__all__ = [
2123
"Constant",
@@ -100,6 +102,22 @@ def __radd__(self, other):
100102
def __rmul__(self, other):
101103
return self.__mul__(other)
102104

105+
def __pow__(self, other):
106+
if(
107+
isinstance(other, theano.compile.SharedVariable) and
108+
other.get_value().squeeze().shape == ()
109+
):
110+
other = tt.squeeze(other)
111+
return Exponentiated(self, other)
112+
elif isinstance(other, Number):
113+
return Exponentiated(self, other)
114+
elif np.asarray(other).squeeze().shape == ():
115+
other = np.squeeze(other)
116+
return Exponentiated(self, other)
117+
118+
raise ValueError("A covariance function can only be exponentiated by a scalar value")
119+
120+
103121
def __array_wrap__(self, result):
104122
"""
105123
Required to allow radd/rmul by numpy arrays.
@@ -172,6 +190,19 @@ def __call__(self, X, Xs=None, diag=False):
172190
return reduce(mul, self.merge_factors(X, Xs, diag))
173191

174192

193+
class Exponentiated(Covariance):
194+
def __init__(self, kernel, power):
195+
self.kernel = kernel
196+
self.power = power
197+
super().__init__(
198+
input_dim=self.kernel.input_dim,
199+
active_dims=self.kernel.active_dims
200+
)
201+
202+
def __call__(self, X, Xs=None, diag=False):
203+
return self.kernel(X, Xs, diag=diag) ** self.power
204+
205+
175206
class Kron(Covariance):
176207
r"""Form a covariance object that is the kronecker product of other covariances.
177208

pymc3/tests/test_gp.py

+55
Original file line numberDiff line numberDiff line change
@@ -237,6 +237,61 @@ def test_multiops(self):
237237
npt.assert_allclose(np.diag(K1), K2d, atol=1e-5)
238238
npt.assert_allclose(np.diag(K2), K1d, atol=1e-5)
239239

240+
class TestCovExponentiation:
241+
def test_symexp_cov(self):
242+
X = np.linspace(0, 1, 10)[:, None]
243+
with pm.Model() as model:
244+
cov1 = pm.gp.cov.ExpQuad(1, 0.1)
245+
cov = cov1 ** 2
246+
K = theano.function([], cov(X))()
247+
npt.assert_allclose(K[0, 1], 0.53940 ** 2, atol=1e-3)
248+
# check diagonal
249+
Kd = theano.function([], cov(X, diag=True))()
250+
npt.assert_allclose(np.diag(K), Kd, atol=1e-5)
251+
252+
def test_covexp_numpy(self):
253+
X = np.linspace(0, 1, 10)[:, None]
254+
with pm.Model() as model:
255+
a = np.array([[2]])
256+
cov = pm.gp.cov.ExpQuad(1, 0.1) ** a
257+
K = theano.function([], cov(X))()
258+
npt.assert_allclose(K[0, 1], 0.53940 ** 2, atol=1e-3)
259+
# check diagonal
260+
Kd = theano.function([], cov(X, diag=True))()
261+
npt.assert_allclose(np.diag(K), Kd, atol=1e-5)
262+
263+
def test_covexp_theano(self):
264+
X = np.linspace(0, 1, 10)[:, None]
265+
with pm.Model() as model:
266+
a = tt.alloc(2.0, 1, 1)
267+
cov = pm.gp.cov.ExpQuad(1, 0.1) ** a
268+
K = theano.function([], cov(X))()
269+
npt.assert_allclose(K[0, 1], 0.53940 ** 2, atol=1e-3)
270+
# check diagonal
271+
Kd = theano.function([], cov(X, diag=True))()
272+
npt.assert_allclose(np.diag(K), Kd, atol=1e-5)
273+
274+
def test_covexp_shared(self):
275+
X = np.linspace(0, 1, 10)[:, None]
276+
with pm.Model() as model:
277+
a = theano.shared(2.0)
278+
cov = pm.gp.cov.ExpQuad(1, 0.1) ** a
279+
K = theano.function([], cov(X))()
280+
npt.assert_allclose(K[0, 1], 0.53940 ** 2, atol=1e-3)
281+
# check diagonal
282+
Kd = theano.function([], cov(X, diag=True))()
283+
npt.assert_allclose(np.diag(K), Kd, atol=1e-5)
284+
285+
def test_invalid_covexp(self):
286+
X = np.linspace(0, 1, 10)[:, None]
287+
with pytest.raises(
288+
ValueError,
289+
match=r"can only be exponentiated by a scalar value"
290+
):
291+
with pm.Model() as model:
292+
a = np.array([[1.0, 2.0]])
293+
cov = pm.gp.cov.ExpQuad(1, 0.1) ** a
294+
240295

241296
class TestCovKron:
242297
def test_symprod_cov(self):

0 commit comments

Comments
 (0)