Package csb :: Package statistics :: Package samplers :: Package mc :: Module singlechain
[frames] | no frames]

Source Code for Module csb.statistics.samplers.mc.singlechain

  1  """ 
  2  Various Monte Carlo equilibrium sampling algorithms, which simulate only one Markov chain. 
  3   
  4  Here is how to sample from a PDF using the L{HMCSampler} class. In the following 
  5  snippet we draw 5000 samples from a 1D normal distribution and plot them: 
  6   
  7   
  8      >>> import numpy 
  9      >>> from csb.io.plots import Chart 
 10      >>> from csb.statistics.pdf import Normal 
 11      >>> from csb.statistics.samplers import State 
 12      >>> from csb.statistics.samplers.mc.singlechain import HMCSampler 
 13       
 14      >>> initial_state = State(numpy.array([1.])) 
 15      >>> grad = lambda q, t: q 
 16      >>> timestep = 1.5 
 17      >>> nsteps = 30 
 18      >>> nsamples = 5000 
 19       
 20      >>> sampler = HMCSampler(Normal(), initial_state, grad, timestep, nsteps) 
 21       
 22      >>> states = [] 
 23      >>> for i in range(nsamples): 
 24              sampler.sample() 
 25              states.append(sampler.state) 
 26       
 27      >>> print('acceptance rate:', sampler.acceptance_rate) 
 28      0.8 
 29       
 30      >>> states = [state.position[0]for state in states] 
 31      >>> chart = Chart() 
 32      >>> chart.plot.hist([numpy.random.normal(size=5000), states], bins=20, normed=True) 
 33      >>> chart.plot.legend(['numpy.random.normal', 'HMC']) 
 34      >>> chart.show() 
 35   
 36   
 37  First, several things which are being needed are imported. 
 38  As every sampler in this module implements a Markov Chain, an initial state has to be 
 39  chosen. In the following lines, several parameters are set: 
 40      - the gradient of the negative log-probability of the PDF under consideration 
 41      - the integration timestep 
 42      - the number of integration steps to be performed in each iteration, that is, the HMC 
 43        trajectory length 
 44      - the number of samples to be drawn 
 45   
 46  The empty list states is initialized. It will serve to store the samples drawn. 
 47  In the loop, C{sampler.sample()} is repeatedly called. After each call of C{sampler.sample()}, 
 48  the current state of the Markov Chain is stored in sampler.state and this state is appended 
 49  to the sample storage list. 
 50   
 51  Then the acceptance rate is printed, the numeric values are being extracted from the 
 52  L{State} objects in states, a histogram is created and finally plotted. 
 53  """ 
 54   
 55  import numpy 
 56  import csb.numeric 
 57   
 58  from csb.statistics.samplers import State 
 59  from csb.statistics.samplers.mc import AbstractSingleChainMC 
 60  from csb.statistics.samplers.mc import SimpleProposalCommunicator 
 61  from csb.statistics.samplers.mc.propagators import MDPropagator 
 62  from csb.numeric.integrators import FastLeapFrog 
63 64 65 -class HMCSampler(AbstractSingleChainMC):
66 """ 67 Hamilton Monte Carlo (HMC, also called Hybrid Monte Carlo by the inventors, 68 Duane, Kennedy, Pendleton, Duncan 1987). 69 70 @param pdf: Probability density function to be sampled from 71 @type pdf: L{csb.statistics.pdf.AbstractDensity} 72 73 @param state: Inital state 74 @type state: L{State} 75 76 @param gradient: Gradient of the negative log-probability 77 @type gradient: L{AbstractGradient} 78 79 @param timestep: Timestep used for integration 80 @type timestep: float 81 82 @param nsteps: Number of integration steps to be performed in 83 each iteration 84 @type nsteps: int 85 86 @param integrator: Subclass of L{AbstractIntegrator} to be used for 87 integrating Hamiltionian equations of motion 88 @type integrator: type 89 """ 90
91 - def __init__(self, pdf, state, gradient, timestep, nsteps, integrator=FastLeapFrog):
92 93 super(HMCSampler, self).__init__(pdf, state) 94 95 self._timestep = None 96 self.timestep = timestep 97 self._nsteps = None 98 self.nsteps = nsteps 99 self._integrator = integrator 100 self._gradient = gradient
101
102 - def _propose(self):
103 104 gen = MDPropagator(self._gradient, self._timestep, self._integrator) 105 momenta = numpy.random.normal(size=self.state.position.shape) 106 self.state = State(self.state.position, momenta) 107 proposal = gen.generate(self.state, self._nsteps).final 108 109 return SimpleProposalCommunicator(proposal)
110
111 - def _calc_pacc(self, proposal_communicator):
112 113 proposal = proposal_communicator.proposal 114 E = lambda x:-self._pdf.log_prob(x) 115 116 pacc = csb.numeric.exp(-0.5 * sum(proposal.momentum ** 2) - E(proposal.position) 117 + 0.5 * sum(self.state.momentum ** 2) + E(self.state.position)) 118 119 return pacc
120 121 @property
122 - def timestep(self):
123 return self._timestep
124 125 @timestep.setter
126 - def timestep(self, value):
127 self._timestep = float(value)
128 129 @property
130 - def nsteps(self):
131 return self._nsteps
132 133 @nsteps.setter
134 - def nsteps(self, value):
135 self._nsteps = int(value)
136
137 -class RWMCSampler(AbstractSingleChainMC):
138 """ 139 Random Walk Metropolis Monte Carlo implementation 140 (Metropolis, Rosenbluth, Teller, Teller 1953; Hastings, 1970). 141 142 @param pdf: Probability density function to be sampled from 143 @type pdf: L{csb.statistics.pdf.AbstractDensity} 144 145 @param state: Inital state 146 @type state: L{State} 147 148 @param stepsize: Serves to set the step size in 149 proposal_density, e.g. for automatic acceptance 150 rate adaption 151 @type stepsize: float 152 153 @param proposal_density: The proposal density as a function f(x, s) 154 of the current state x and the stepsize s. 155 By default, the proposal density is uniform, 156 centered around x, and has width s. 157 @type proposal_density: callable 158 """ 159
160 - def __init__(self, pdf, state, stepsize=1., proposal_density=None):
161 162 super(RWMCSampler, self).__init__(pdf, state) 163 self._stepsize = None 164 self.stepsize = stepsize 165 if proposal_density == None: 166 self._proposal_density = lambda x, s: x.position + s * numpy.random.uniform(size=x.position.shape, low=-1., high=1.) 167 else: 168 self._proposal_density = proposal_density
169
170 - def _propose(self):
171 172 return SimpleProposalCommunicator(State(self._proposal_density(self._state, self.stepsize)))
173
174 - def _calc_pacc(self, proposal_communicator):
175 176 proposal = proposal_communicator.proposal 177 E = lambda x:-self._pdf.log_prob(x) 178 179 pacc = csb.numeric.exp(-(E(proposal.position) - E(self.state.position))) 180 return pacc
181 182 @property
183 - def stepsize(self):
184 return self._stepsize
185 186 @stepsize.setter
187 - def stepsize(self, value):
188 self._stepsize = float(value)
189