gstlal-inspiral  0.4.2
 All Classes Namespaces Files Functions Variables Pages
ptsampler.py
1 #!/usr/bin/env python
2 # -*- coding: utf-8 -*-
3 
4 from __future__ import (division, print_function, absolute_import,
5  unicode_literals)
6 
7 __all__ = ["PTSampler"]
8 
9 try:
10  import acor
11  acor = acor
12 except ImportError:
13  acor = None
14 
15 from .sampler import Sampler
16 from .ensemble import EnsembleSampler
17 import multiprocessing as multi
18 import numpy as np
19 import numpy.random as nr
20 
21 
22 class PTPost(object):
23  """
24  Wrapper for posterior used with the :class:`PTSampler` in emcee.
25 
26  """
27 
28  def __init__(self, logl, logp, beta):
29  """
30  :param logl:
31  Function returning natural log of the likelihood.
32 
33  :param logp:
34  Function returning natural log of the prior.
35 
36  :param beta:
37  Inverse temperature of this chain: ``lnpost = beta*logl + logp``.
38 
39  """
40 
41  self._logl = logl
42  self._logp = logp
43  self._beta = beta
44 
45  def __call__(self, x):
46  """
47  Returns ``lnpost(x)``, ``lnlike(x)`` (the second value will be
48  treated as a blob by emcee), where
49 
50  .. math::
51 
52  \ln \pi(x) \equiv \beta \ln l(x) + \ln p(x)
53 
54  :param x:
55  The position in parameter space.
56 
57  """
58  lp = self._logp(x)
59 
60  # If outside prior bounds, return 0.
61  if lp == float('-inf'):
62  return lp, lp
63 
64  ll = self._logl(x)
65 
66  return self._beta * ll + lp, ll
67 
68 
69 class PTSampler(Sampler):
70  """
71  A parallel-tempered ensemble sampler, using :class:`EnsembleSampler`
72  for sampling within each parallel chain.
73 
74  :param ntemps:
75  The number of temperatures.
76 
77  :param nwalkers:
78  The number of ensemble walkers at each temperature.
79 
80  :param dim:
81  The dimension of parameter space.
82 
83  :param logl:
84  The log-likelihood function.
85 
86  :param logp:
87  The log-prior function.
88 
89  :param threads: (optional)
90  The number of parallel threads to use in sampling.
91 
92  :param pool: (optional)
93  Alternative to ``threads``. Any object that implements a
94  ``map`` method compatible with the built-in ``map`` will do
95  here. For example, :class:`multi.Pool` will do.
96 
97  :param betas: (optional)
98  Array giving the inverse temperatures, :math:`\\beta=1/T`, used in the
99  ladder. The default is for an exponential ladder, with beta
100  decreasing by a factor of :math:`1/\\sqrt{2}` each rung.
101 
102  """
103  def __init__(self, ntemps, nwalkers, dim, logl, logp, threads=1,
104  pool=None, betas=None):
105  self.logl = logl
106  self.logp = logp
107 
108  self.ntemps = ntemps
109  self.nwalkers = nwalkers
110  self.dim = dim
111 
112  self._chain = None
113  self._lnprob = None
114  self._lnlikelihood = None
115 
116  if betas is None:
117  self._betas = self.exponential_beta_ladder(ntemps)
118  else:
119  self._betas = betas
120 
121  self.nswap = np.zeros(ntemps, dtype=np.float)
122  self.nswap_accepted = np.zeros(ntemps, dtype=np.float)
123 
124  self.pool = pool
125  if threads > 1 and pool is None:
126  self.pool = multi.Pool(threads)
127 
128  self.samplers = [EnsembleSampler(nwalkers, dim,
129  PTPost(logl, logp, b),
130  pool=self.pool)
131  for b in self.betas]
132 
133  def exponential_beta_ladder(self, ntemps):
134  """
135  Exponential ladder in :math:`1/T`, with :math:`T` increasing by
136  :math:`\\sqrt{2}` each step, with ``ntemps`` in total.
137 
138  """
139  return np.exp(np.linspace(0, -(ntemps - 1) * 0.5 * np.log(2), ntemps))
140 
141  def reset(self):
142  """
143  Clear the ``chain``, ``lnprobability``, ``lnlikelihood``,
144  ``acceptance_fraction``, ``tswap_acceptance_fraction`` stored
145  properties.
146 
147  """
148  for s in self.samplers:
149  s.reset()
150 
151  self.nswap = np.zeros(self.ntemps, dtype=np.float)
152  self.nswap_accepted = np.zeros(self.ntemps, dtype=np.float)
153 
154  self._chain = None
155  self._lnprob = None
156  self._lnlikelihood = None
157 
158  def sample(self, p0, lnprob0=None, lnlike0=None, iterations=1,
159  thin=1, storechain=True):
160  """
161  Advance the chains ``iterations`` steps as a generator.
162 
163  :param p0:
164  The initial positions of the walkers. Shape should be
165  ``(ntemps, nwalkers, dim)``.
166 
167  :param lnprob0: (optional)
168  The initial posterior values for the ensembles. Shape
169  ``(ntemps, nwalkers)``.
170 
171  :param lnlike0: (optional)
172  The initial likelihood values for the ensembles. Shape
173  ``(ntemps, nwalkers)``.
174 
175  :param iterations: (optional)
176  The number of iterations to preform.
177 
178  :param thin: (optional)
179  The number of iterations to perform between saving the
180  state to the internal chain.
181 
182  :param storechain: (optional)
183  If ``True`` store the iterations in the ``chain``
184  property.
185 
186  At each iteration, this generator yields
187 
188  * ``p``, the current position of the walkers.
189 
190  * ``lnprob`` the current posterior values for the walkers.
191 
192  * ``lnlike`` the current likelihood values for the walkers.
193 
194  """
195  p = np.copy(np.array(p0))
196 
197  # If we have no lnprob or logls compute them
198  if lnprob0 is None or lnlike0 is None:
199  lnprob0 = np.zeros((self.ntemps, self.nwalkers))
200  lnlike0 = np.zeros((self.ntemps, self.nwalkers))
201  for i in range(self.ntemps):
202  fn = PTPost(self.logl, self.logp, self.betas[i])
203  if self.pool is None:
204  results = list(map(fn, p[i, :, :]))
205  else:
206  results = list(self.pool.map(fn, p[i, :, :]))
207 
208  lnprob0[i, :] = np.array([r[0] for r in results])
209  lnlike0[i, :] = np.array([r[1] for r in results])
210 
211  lnprob = lnprob0
212  logl = lnlike0
213 
214  # Expand the chain in advance of the iterations
215  if storechain:
216  nsave = iterations / thin
217  if self._chain is None:
218  isave = 0
219  self._chain = np.zeros((self.ntemps, self.nwalkers, nsave,
220  self.dim))
221  self._lnprob = np.zeros((self.ntemps, self.nwalkers, nsave))
222  self._lnlikelihood = np.zeros((self.ntemps, self.nwalkers,
223  nsave))
224  else:
225  isave = self._chain.shape[2]
226  self._chain = np.concatenate((self._chain,
227  np.zeros((self.ntemps,
228  self.nwalkers,
229  nsave, self.dim))),
230  axis=2)
231  self._lnprob = np.concatenate((self._lnprob,
232  np.zeros((self.ntemps,
233  self.nwalkers,
234  nsave))),
235  axis=2)
236  self._lnlikelihood = np.concatenate((self._lnlikelihood,
237  np.zeros((self.ntemps,
238  self.nwalkers,
239  nsave))),
240  axis=2)
241 
242  for i in range(iterations):
243  for j, s in enumerate(self.samplers):
244  for psamp, lnprobsamp, rstatesamp, loglsamp in s.sample(
245  p[j, ...],
246  lnprob0=lnprob[j, ...],
247  blobs0=logl[j, ...], storechain=False):
248  p[j, ...] = psamp
249  lnprob[j, ...] = lnprobsamp
250  logl[j, ...] = np.array(loglsamp)
251 
252  p, lnprob, logl = self._temperature_swaps(p, lnprob, logl)
253 
254  if (i + 1) % thin == 0:
255  if storechain:
256  self._chain[:, :, isave, :] = p
257  self._lnprob[:, :, isave, ] = lnprob
258  self._lnlikelihood[:, :, isave] = logl
259  isave += 1
260 
261  yield p, lnprob, logl
262 
263  def _temperature_swaps(self, p, lnprob, logl):
264  """
265  Perform parallel-tempering temperature swaps on the state
266  in ``p`` with associated ``lnprob`` and ``logl``.
267 
268  """
269  ntemps = self.ntemps
270 
271  for i in range(ntemps - 1, 0, -1):
272  bi = self.betas[i]
273  bi1 = self.betas[i - 1]
274 
275  dbeta = bi1 - bi
276 
277  iperm = nr.permutation(self.nwalkers)
278  i1perm = nr.permutation(self.nwalkers)
279 
280  raccept = np.log(nr.uniform(size=self.nwalkers))
281  paccept = dbeta * (logl[i, iperm] - logl[i - 1, i1perm])
282 
283  self.nswap[i] += self.nwalkers
284  self.nswap[i - 1] += self.nwalkers
285 
286  asel = (paccept > raccept)
287  nacc = np.sum(asel)
288 
289  self.nswap_accepted[i] += nacc
290  self.nswap_accepted[i - 1] += nacc
291 
292  ptemp = np.copy(p[i, iperm[asel], :])
293  ltemp = np.copy(logl[i, iperm[asel]])
294  prtemp = np.copy(lnprob[i, iperm[asel]])
295 
296  p[i, iperm[asel], :] = p[i - 1, i1perm[asel], :]
297  logl[i, iperm[asel]] = logl[i - 1, i1perm[asel]]
298  lnprob[i, iperm[asel]] = lnprob[i - 1, i1perm[asel]] \
299  - dbeta * logl[i - 1, i1perm[asel]]
300 
301  p[i - 1, i1perm[asel], :] = ptemp
302  logl[i - 1, i1perm[asel]] = ltemp
303  lnprob[i - 1, i1perm[asel]] = prtemp + dbeta * ltemp
304 
305  return p, lnprob, logl
306 
307  def thermodynamic_integration_log_evidence(self, logls=None, fburnin=0.1):
308  """
309  Thermodynamic integration estimate of the evidence.
310 
311  :param logls: (optional) The log-likelihoods to use for
312  computing the thermodynamic evidence. If ``None`` (the
313  default), use the stored log-likelihoods in the sampler.
314  Should be of shape ``(Ntemps, Nwalkers, Nsamples)``.
315 
316  :param fburnin: (optional)
317  The fraction of the chain to discard as burnin samples; only the
318  final ``1-fburnin`` fraction of the samples will be used to
319  compute the evidence; the default is ``fburnin = 0.1``.
320 
321  :return ``(lnZ, dlnZ)``: Returns an estimate of the
322  log-evidence and the error associated with the finite
323  number of temperatures at which the posterior has been
324  sampled.
325 
326  The evidence is the integral of the un-normalized posterior
327  over all of parameter space:
328 
329  .. math::
330 
331  Z \\equiv \\int d\\theta \\, l(\\theta) p(\\theta)
332 
333  Thermodymanic integration is a technique for estimating the
334  evidence integral using information from the chains at various
335  temperatures. Let
336 
337  .. math::
338 
339  Z(\\beta) = \\int d\\theta \\, l^\\beta(\\theta) p(\\theta)
340 
341  Then
342 
343  .. math::
344 
345  \\frac{d \\ln Z}{d \\beta}
346  = \\frac{1}{Z(\\beta)} \\int d\\theta l^\\beta p \\ln l
347  = \\left \\langle \\ln l \\right \\rangle_\\beta
348 
349  so
350 
351  .. math::
352 
353  \\ln Z(\\beta = 1)
354  = \\int_0^1 d\\beta \\left \\langle \\ln l \\right\\rangle_\\beta
355 
356  By computing the average of the log-likelihood at the
357  difference temperatures, the sampler can approximate the above
358  integral.
359  """
360 
361  if logls is None:
363  logls=self.lnlikelihood, fburnin=fburnin)
364  else:
365  betas = np.concatenate((self.betas, np.array([0])))
366  betas2 = np.concatenate((self.betas[::2], np.array([0])))
367 
368  istart = int(logls.shape[2] * fburnin + 0.5)
369 
370  mean_logls = np.mean(np.mean(logls, axis=1)[:, istart:], axis=1)
371  mean_logls2 = mean_logls[::2]
372 
373  lnZ = -np.dot(mean_logls, np.diff(betas))
374  lnZ2 = -np.dot(mean_logls2, np.diff(betas2))
375 
376  return lnZ, np.abs(lnZ - lnZ2)
377 
378  @property
379  def betas(self):
380  """
381  Returns the sequence of inverse temperatures in the ladder.
382 
383  """
384  return self._betas
385 
386  @property
387  def chain(self):
388  """
389  Returns the stored chain of samples; shape ``(Ntemps,
390  Nwalkers, Nsteps, Ndim)``.
391 
392  """
393  return self._chain
394 
395  @property
396  def lnprobability(self):
397  """
398  Matrix of lnprobability values; shape ``(Ntemps, Nwalkers, Nsteps)``.
399 
400  """
401  return self._lnprob
402 
403  @property
404  def lnlikelihood(self):
405  """
406  Matrix of ln-likelihood values; shape ``(Ntemps, Nwalkers, Nsteps)``.
407 
408  """
409  return self._lnlikelihood
410 
411  @property
413  """
414  Returns an array of accepted temperature swap fractions for
415  each temperature; shape ``(ntemps, )``.
416 
417  """
418  return self.nswap_accepted / self.nswap
419 
420  @property
422  """
423  Matrix of shape ``(Ntemps, Nwalkers)`` detailing the
424  acceptance fraction for each walker.
425 
426  """
427  return np.array([s.acceptance_fraction for s in self.samplers])
428 
429  @property
430  def acor(self):
431  """
432  Returns a matrix of autocorrelation lengths for each
433  parameter in each temperature of shape ``(Ntemps, Ndim)``.
434 
435  """
436  if acor is None:
437  raise ImportError('acor')
438  else:
439  acors = np.zeros((self.ntemps, self.dim))
440 
441  for i in range(self.ntemps):
442  for j in range(self.dim):
443  acors[i, j] = acor.acor(self._chain[i, :, :, j])[0]
444 
445  return acors