Models

Model object

class ddm.model.Model(drift=DriftConstant(drift=0), noise=NoiseConstant(noise=1), bound=BoundConstant(B=1), IC=ICPointSourceCenter(), overlay=OverlayNone(), name='', dx=0.005, dt=0.005, T_dur=2.0, fitresult=None)[source]

Bases: object

A full simulation of a single DDM-style model.

Each model simulation depends on five key components:

  • A description of how drift rate (drift) changes throughout the simulation.
  • A description of how variability (noise) changes throughout the simulation.
  • A description of how the boundary changes throughout the simulation.
  • Starting conditions for the model
  • Specific details of a task which cause dynamic changes in the model (e.g. a stimulus intensity change)

This class manages these, and also provides the affiliated services, such as analytical or numerical simulations of the resulting reaction time distribution.

IC(conditions)[source]

The initial distribution at t=0.

Returns a length N ndarray (where N is the size of x_domain()) which should sum to 1.

can_solve_cn(conditions={})[source]

Check whether this model is compatible with Crank-Nicolson solver.

All bound functions which do not depend on time are compatible.

can_solve_explicit(conditions={})[source]

Check explicit method stability criterion

flux(x, t, conditions)[source]

The flux across the boundary at position x at time t.

get_dependence(name)[source]

Return the dependence object given by the string name.

get_fit_result()[source]

Returns a FitResult object describing how the model was fit.

Returns the FitResult object describing the last time this model was fit to data, including the loss function, fitting method, and the loss function value. If the model was never fit to data, this will return FitResultEmpty.

get_model_parameter_names()[source]

Get an ordered list of the names of all parameters in the model.

Returns the name of each model parameter. The ordering is arbitrary, but is uaranteed to be in the same order as get_model_parameters() and set_model_parameters(). If multiple parameters refer to the same “Fittable” object, then that object will only be listed once, however the names of the parameters will be separated by a “/” character.

get_model_parameters()[source]

Get an ordered list of all model parameters.

Returns a list of each model parameter which can be varied during a fitting procedure. The ordering is arbitrary but is guaranteed to be in the same order as get_model_parameter_names() and set_model_parameters(). If multiple parameters refer to the same “Fittable” object, then that object will only be listed once.

get_model_type()[source]

Return a dictionary which fully specifies the class of the five key model components.

has_analytical_solution()[source]

Is it possible to find an analytic solution for this model?

parameters()[source]

Return all parameters in the model

This will return a dictionary of dictionaries. The keys of this dictionary will be “drift”, “noise”, “bound”, “IC”, and “overlay”. The values will be dictionaries, each containing the parameters used by these. Note that this includes both fixed parameters and Fittable parameters. If a parameter is fittable, it will return either a Fittable object or a Fitted object in place of the parameter, depending on whether or not the model has been fit to data yet.

set_model_parameters(params)[source]

Set the parameters of the model from an ordered list.

Takes as an argument a list of parameters in the same order as those from get_model_parameters(). Sets the associated parameters as a “Fitted” object. If multiple parameters refer to the same “Fittable” object, then that object will only be listed once.

simulate_trial(conditions={}, cutoff=True, rk4=True, seed=0)[source]

Simulate the decision variable for one trial.

Given conditions conditions, this function will simulate the decision variable for a single trial. It will cut off the simulation when the decision variable crosses the boundary unless cutoff is set to False. By default, Runge-Kutta is used to simulate the trial, however if rk4 is set to False, the less efficient Euler’s method is used instead. This returns a trajectory of the simulated trial over time as a numpy array.

Note that this will return the same trajectory on each run unless the random seed seed is varied.

Also note that you shouldn’t normally need to use this function. To simulate an entire probability distributions, call Model.solve() and the results of the simulation will be in the returned Solution object. This is only useful for finding individual trajectories instead of the probability distribution as a whole.

simulated_solution(conditions={}, size=1000, rk4=True, seed=0)[source]

Simulate individual trials to obtain a distribution.

Given conditions conditions and the number size of trials to simulate, this will run the function “simulate_trial” size times, and use the result to find a histogram analogous to solve. Returns a Sample object.

Note that in practice you should never need to use this function. This function uses an outdated method to simulate the model and should be used for comparison perposes only. To produce a probability density function of boundary crosses, use Model.solve(). To sample from the probability distribution (e.g. for finding confidence intervals for limited amounts of data), call Model.solve() and then use the Solution.resample() function of the resulting Solution.

solve(conditions={}, return_evolution=False, force_python=False)[source]

Solve the model using an analytic solution if possible, and a numeric solution if not.

First, it tries to use Crank-Nicolson as the solver, and then backward Euler. See documentation of Model.solve_numerical() for more information.

The return_evolution argument should be set to True if you need to use the Solution.get_evolution() function from the returned Solution.

Return a Solution object describing the joint PDF distribution of reaction times.

solve_analytical(conditions={}, force_python=False)[source]

Solve the model with an analytic solution, if possible.

Analytic solutions are only possible in a select number of special cases; in particular, it works for simple DDM and for linearly collapsing bounds and arbitrary single-point initial conditions. (See Anderson (1960) for implementation details.) For most reasonably complex models, the method will fail. Check whether a solution is possible with has_analytic_solution().

If successful, this returns a Solution object describing the joint PDF. If unsuccessful, this will raise an exception.

solve_numerical(method='cn', conditions={}, return_evolution=False, force_python=False)[source]

Solve the DDM model numerically.

Use method to solve the DDM. method can either be “explicit”, “implicit”, or “cn” (for Crank-Nicolson). This is the core DDM solver of this library.

Crank-Nicolson is the default and works for any model with constant bounds.

Implicit is the fallback method. It should work well in most cases and is generally stable.

Normally, the explicit method should not be used. Also note the stability criteria for explicit method is:

noise^2/2 * dt/dx^2 < 1/2

It returns a Solution object describing the joint PDF. This method should not fail for any model type.

return_evolution (default=False) governs whether or not the function returns the full evolution of the pdf as part of the Solution object. This only works with methods “explicit” or “implicit”, not with “cn”.

force_python makes PyDDM use the solver written in Python instead of the optimized solver written in C.

solve_numerical_c(conditions={})[source]

Solve the DDM model using the implicit method with C extensions.

This function should give near identical results to solve_numerical_implicit. However, it uses compiled C code instead of Python code to do so, which should make it much (10-100x) faster.

This does not current work with non-Gaussian diffusion matrices (a currently undocumented feature).

solve_numerical_cn(conditions={})[source]

Solve the DDM model numerically using Crank-Nicolson.

This uses the Crank Nicolson method to solve the DDM at each timepoint. Results are then compiled together. This is the core DDM solver of this library.

It returns a Solution object describing the joint PDF.

solve_numerical_explicit(conditions={}, **kwargs)[source]

Solve the model using the explicit method (Forward Euler).

See documentation for the solve_numerical method.

solve_numerical_implicit(conditions={}, **kwargs)[source]

Solve the model using the implicit method (Backward Euler).

See documentation for the solve_numerical method.

t_domain()[source]

A list of all of the timepoints over which the joint PDF will be defined (increments of dt from 0 to T_dur).

x_domain(conditions, t=None)[source]

A list which spans from the lower boundary to the upper boundary by increments of dx.

Fittable object

class ddm.model.Fittable(**kwargs)[source]

Bases: float

For parameters that should be adjusted when fitting a model to data.

Each Fittable object does not need any parameters, however several parameters may improve the ability to fit the model. In particular, maxval and minval ensure we do not choose an invalid parameter value. default is the value to start with when fitting; if it is not given, it will be selected at random.

default()[source]

Choose a default value.

This chooses a value for the Fittable object randomly abiding by any constraints. Note that calling this function multiple times will give different results.

Fitted object

class ddm.model.Fitted(val, **kwargs)[source]

Bases: ddm.model.Fittable

Parameters which used to be Fittable but now hold a value.

This extends float so that the parameters for the Fittable can be saved in the final model, even though there are no more Fittable objects in the final model.

Sample object (for RT data)

Solution object (for theoretical RT histograms)