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

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

  1  """ 
  2  Provides various deterministic and stochastic propagators. 
  3  """ 
  4   
  5  import numpy 
  6   
  7  from abc import ABCMeta, abstractmethod 
  8   
  9  from csb.statistics.samplers import State 
 10  from csb.statistics.samplers.mc import TrajectoryBuilder 
 11  from csb.numeric.integrators import FastLeapFrog, VelocityVerlet 
12 -class AbstractPropagator(object):
13 """ 14 Abstract propagator class. Subclasses serve to propagate 15 an inital state by some dynamics to a final state. 16 """ 17 18 __metaclass__ = ABCMeta 19 20 @abstractmethod
21 - def generate(self, init_state, length, return_trajectory=False):
22 """ 23 Generate a trajectory, starting from an initial state with a certain length. 24 25 @param init_state: Initial state from which to propagate 26 @type init_state: L{State} 27 28 @param length: Length of the trajectory (in integration steps or stochastic moves) 29 @type length: int 30 31 @param return_trajectory: Return complete L{Trajectory} instead of the initial 32 and final states only (L{PropagationResult}) 33 @type return_trajectory: boolean 34 35 @rtype: L{AbstractPropagationResult} 36 """ 37 pass
38
39 -class MDPropagator(AbstractPropagator):
40 """ 41 Molecular Dynamics propagator. Generates a trajectory 42 by integration of Hamiltionian equations of motion. 43 44 @param gradient: Gradient of potential energy. Guides the dynamics. 45 @type gradient: L{AbstractGradient} 46 47 @param timestep: Timestep to be used for integration 48 @type timestep: float 49 50 @param integrator: Subclass of L{AbstractIntegrator} to be used to integrate 51 Hamiltonian equations of motion 52 @type integrator: type 53 """ 54
55 - def __init__(self, gradient, timestep, integrator=FastLeapFrog):
56 57 super(MDPropagator, self).__init__() 58 59 self._gradient = None 60 self.gradient = gradient 61 62 self._timestep = None 63 self.timestep = timestep 64 65 self._integrator = integrator
66 67 @property
68 - def gradient(self):
69 return self._gradient
70 @gradient.setter
71 - def gradient(self, value):
72 self._gradient = value
73 74 @property
75 - def timestep(self):
76 return self._timestep
77 @timestep.setter
78 - def timestep(self, value):
79 self._timestep = float(value)
80
81 - def generate(self, init_state, length, return_trajectory=False):
82 83 integrator = self._integrator(self.timestep, self.gradient) 84 85 result = integrator.integrate(init_state, length, return_trajectory) 86 return result
87
88 -class ThermostattedMDPropagator(MDPropagator):
89 """ 90 Thermostatted Molecular Dynamics propagator. Employs the Andersen thermostat 91 which simulates collision with particles of a heat bath at a given temperature. 92 93 @param gradient: Gradient of potential energy. Guides the dynamics. 94 @type gradient: L{AbstractGradient} 95 96 @param timestep: Timestep to be used for integration 97 @type timestep: float 98 99 @param temperature: Time-dependent temperature 100 @type temperature: Real-valued function 101 102 @param collision_probability: collision probability within duration of one timestep 103 @type collision_probability: float 104 105 @param update_interval: Interval with which momenta are redrawn 106 @type update_interval: int 107 108 @param integrator: Subclass of L{AbstractIntegrator} to be used to perform 109 integration steps between momentum updates 110 @type integrator: type 111 """ 112
113 - def __init__(self, gradient, timestep, temperature=lambda t: 1., collision_probability=0.1, 114 update_interval=1, integrator=VelocityVerlet):
115 116 super(ThermostattedMDPropagator, self).__init__(gradient, timestep, integrator) 117 118 self._collision_probability = collision_probability 119 self._update_interval = update_interval 120 self._temperature = temperature
121
122 - def _update(self, momentum, T, collision_probability):
123 """ 124 Simulate collision with heat bath particles. 125 126 @param momentum: Momentum 127 @type momentum: one-dimensional numpy array of numbers 128 129 @param T: Temperature of the heat bath 130 @type T: float 131 132 @param collision_probability: collision probability within duration of one timestep 133 @type collision_probability: float 134 135 @rtype: tuple (updated momentum, heat induced by the update) 136 """ 137 138 heat = 0. 139 update_list = [] 140 for k in range(len(momentum)): 141 if numpy.random.random() < collision_probability: 142 update_list.append(k) 143 if len(update_list) > 0: 144 ke_old = 0.5 * sum(momentum ** 2) 145 for k in range(len(momentum)): 146 if k in update_list: momentum[k] = numpy.random.normal(scale=numpy.sqrt(T)) 147 heat = (0.5 * sum(momentum ** 2) - ke_old) / T 148 149 return momentum, heat
150
151 - def _step(self, i, state, heat, integrator):
152 """ 153 Performs one step consisting of an integration step 154 and possibly a momentum update 155 156 @param i: integration step count 157 @type i: int 158 159 @param state: state to be updated 160 @type state: L{State} 161 162 @param heat: heat produced up to the current integration step 163 @type heat: float 164 165 @param integrator: integration scheme used to evolve the state deterministically 166 @type integrator: L{AbstractIntegrator} 167 """ 168 state = integrator.integrate_once(state, i) 169 170 if i % self._update_interval == 0: 171 state.momentum, stepheat = self._update(state.momentum, 172 self._temperature(i * self.timestep), 173 self._collision_probability) 174 175 heat += stepheat 176 177 return state, heat
178
179 - def generate(self, init_state, length, return_trajectory=False):
180 181 integrator = self._integrator(self.timestep, self.gradient) 182 builder = TrajectoryBuilder.create(full=return_trajectory) 183 184 builder.add_initial_state(init_state) 185 186 heat = 0. 187 state = init_state.clone() 188 189 for i in range(length - 1): 190 state, heat = self._step(i, state, heat, integrator) 191 builder.add_intermediate_state(state) 192 193 state, heat = self._step(length - 1, state, heat, integrator) 194 builder.add_final_state(state) 195 196 traj = builder.product 197 traj.heat = heat 198 199 return traj
200
201 -class AbstractMCPropagator(AbstractPropagator):
202 """ 203 Provides the interface for MC trajectory generators. Implementations 204 generate a sequence of states according to some implementation of 205 L{AbstractSingleChainMC}. 206 207 @param pdf: PDF to sample from 208 @type pdf: L{AbstractDensity} 209 210 @param temperature: See documentation of L{AbstractSingleChainMC} 211 @type temperature: float 212 """ 213 214 __metaclass__ = ABCMeta 215
216 - def __init__(self, pdf, temperature=1.):
217 218 self._pdf = pdf 219 self._temperature = temperature 220 self._acceptance_rate = 0.0
221
222 - def generate(self, init_state, length, return_trajectory=True):
223 224 self._init_sampler() 225 self._sampler.state = init_state 226 227 builder = TrajectoryBuilder.create(full=return_trajectory) 228 229 builder.add_initial_state(init_state) 230 231 for i in range(length): 232 self._sampler.sample() 233 builder.add_intermediate_state(self._sampler.state) 234 235 self._acceptance_rate = self._sampler.acceptance_rate 236 237 return builder.product
238 239 @abstractmethod
240 - def _init_sampler(self):
241 """ 242 Initializes the sampler with which to obtain the MC state 243 trajectory. 244 """ 245 246 pass
247 248 @property
249 - def acceptance_rate(self):
250 """ 251 Acceptance rate of the MC sampler that generated the 252 trajectory. 253 """ 254 return self._acceptance_rate
255
256 -class RWMCPropagator(AbstractMCPropagator):
257 """ 258 Draws a number of samples from a PDF using the L{RWMCSampler} and 259 returns them as a L{Trajectory}. 260 261 @param pdf: PDF to sample from 262 @type pdf: L{AbstractDensity} 263 @param stepsize: Serves to set the step size in 264 proposal_density, e.g. for automatic acceptance 265 rate adaption 266 @type stepsize: float 267 @param proposal_density: The proposal density as a function f(x, s) 268 of the current state x and the stepsize s. 269 By default, the proposal density is uniform, 270 centered around x, and has width s. 271 @type proposal_density: callable 272 273 @param temperature: See documentation of L{AbstractSingleChainMC} 274 @type temperature: float 275 """ 276
277 - def __init__(self, pdf, stepsize=1., proposal_density=None, temperature=1.):
278 279 super(RWMCPropagator, self).__init__(pdf, temperature) 280 281 self._stepsize = stepsize 282 self._proposal_density = proposal_density 283 284 self._init_sampler()
285
286 - def _init_sampler(self):
287 288 from csb.statistics.samplers.mc.singlechain import RWMCSampler 289 290 dummy_state = State(numpy.array([0.])) 291 self._sampler = RWMCSampler(self._pdf, dummy_state, self._stepsize, 292 self._proposal_density, self._temperature)
293
294 -class HMCPropagator(AbstractMCPropagator):
295 """ 296 Draws a number of samples from a PDF using the L{HMCSampler} and 297 returns them as a L{Trajectory}. 298 299 @param pdf: PDF to sample from 300 @type pdf: L{AbstractDensity} 301 @param gradient: Gradient of the negative log-probability 302 @type gradient: L{AbstractGradient} 303 304 @param timestep: Timestep used for integration 305 @type timestep: float 306 307 @param nsteps: Number of integration steps to be performed in 308 each iteration 309 @type nsteps: int 310 311 @param integrator: Subclass of L{AbstractIntegrator} to be used for 312 integrating Hamiltionian equations of motion 313 @type integrator: type 314 315 @param temperature: See documentation of L{AbstractSingleChainMC} 316 @type temperature: float 317 """ 318
319 - def __init__(self, pdf, gradient, timestep, nsteps, integrator=FastLeapFrog, temperature=1.):
320 321 super(HMCPropagator, self).__init__(pdf, temperature) 322 323 self._gradient = gradient 324 self._timestep = timestep 325 self._nsteps = nsteps 326 self._integrator = integrator
327
328 - def _init_sampler(self):
329 330 from csb.statistics.samplers.mc.singlechain import HMCSampler 331 332 dummy_state = State(numpy.array([0.])) 333 self._sampler = HMCSampler(self._pdf, dummy_state, self._gradient, 334 self._timestep, self._nsteps, self._integrator, 335 self._temperature)
336