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