OpenCL Simulation

Myokit provides a simulation engine for parallelized multi-cellular simulations, implemented using OpenCL. The engine can handle 1d and 2d rectangular tissue or networks of arbitrarily connected cells. Heterogeneity can be introduced by replacing model constants with fields, using the method set_field.

A two-model fiber and tissue simulation can be run using the FiberTissueSimulation.

Information about the available OpenCL devices can be obtained using the OpenCL class, which also allows the preferred device to be selected. Note that this functionality is also accessible through the command-line myo script (see opencl-select).

class myokit.SimulationOpenCL(model, protocol=None, ncells=256, diffusion=True, precision=32, native_maths=False, rl=False)

Can run multi-cellular simulations based on a model using OpenCL for parallelization.

Takes the following input arguments:

model

The model to simulate with. This model will be cloned when the simulation is created so that no changes to the given model will be made.

protocol

An optional pacing protocol, used to stimulate a number of cells either at the start of a fiber or at the bottom-left of the tissue.

ncells

The number of cells. Use a scalar for 1d simulations or a tuple (nx, ny) for 2d simulations.

diffusion

Can be set to False to disable diffusion currents. This can be useful in combination with set_field() to explore the effects of varying one or more parameters in a single cell model.

precision

Can be set to myokit.SINGLE_PRECISION (default) or myokit.DOUBLE_PRECISION if the used device supports it.

native_maths

On some devices, selected functions (e.g. exp) can be made to run faster (but possibly less accurately) by setting native_maths=True.

rl

Use Rush-Larsen updates instead of forward Euler for any Hodgkin-Huxley gating variables (default=``False``).

The simulation provides the following inputs variables can bind to:

time

The simulation time

pace

The pacing level, this is set if a protocol was passed in. Will be set to 0 if no protocol is provided.

diffusion_current (if enabled)

The current flowing from the cell to its neighbors. This will be positive when the cell is acting as a source, negative when it is acting as a sink.

The input time is set globally: Any variable bound to time will appear in the logs as single, global variable (for example engine.time instead of 1.2.engine.time. Variables bound to pace or diffusion_current are logged per cell. (If diffusion currents are disabled, the input diffusion_current will not be used, and any variables bound to it will be logged or not according to their default value).

To set the number of cells that will be paced, the methods set_paced_cells() and set_paced_cell_list() can be used.

A single labeled variable is required for this simulation to work:

membrane_potential

The variable representing the membrane potential.

Simulations maintain an internal state consisting of

  • the current simulation time

  • the current state

  • the default state

When a simulation is created, the simulation time is set to 0 and both the current and the default state are equal to the state of the given model, copied once for each cell. After each call to run() the time variable and current state are updated, so that each successive call to run continues where the previous simulation left off. A reset() method is provided that will set the time back to 0 and revert the current state to the default state. To change the time or state manually, use set_time() and set_state().

A pre-pacing method pre() is provided that doesn’t affect the simulation time but will update the current and the default state. This allows you to pre-pace, run a simulation, reset to the pre-paced state, run another simulation etc.

To set up a 1d simulation, the argument ncells should be given as a tuple. In this case, any cell i will be assumed to be connected to cells i - 1 and i + 1 (except at the boundaries). Similarly, if ncells is a 2-dimensional tuple (nx, ny) a grid will be assumed so that each (non-boundary) cell is connected to four neighbors. Finally, arbitrary geometries can be used by passing a scalar to ncells and specifying the connections with set_connections().

The diffusion_current in any cell is calculated as:

i = sum[g * (V - V_j)]

Where the sum is taken over all connected cells j (see [1]).

Models used with this simulation need to have independent components: it should be possible to evaluate the model’s equations one component at a time. A model’s suitability can be tested using has_interdependent_components.

Simulations are performed using a forward Euler (FE) method with a fixed time step (see set_step_size()). Using smaller step sizes is computationally demanding, but gives more accurate results. Using too large a step size can also cause a simulation to become unstable, but please note that stability does imply accuracy, and always double-check important results by re-running with a reduced step size.

If the optional parameter rl is set to True, state variables written in a Hodgkin-Huxley form will be updated using a Rush-Larsen (RL) instead of a forward Euler step (see [2]). This provides greater stability (so that the step size can be increased) but not necessarily greater accuracy (see [3]), so that care must be taken when using this method.

[1] Myokit: A simple interface to cardiac cellular electrophysiology. Clerx, Collins, de Lange, Volders (2016) Progress in Biophysics and Molecular Biology.

[2] A practical algorithm for solving dynamic membrane equations. Rush, Larsen (1978) IEEE Transactions on Biomedical Engineering

[3] Cellular cardiac electrophysiology modelling with Chaste and CellML Cooper, Spiteri, Mirams (2015) Frontiers in Physiology

calculate_conductance(r, sx, chi, dx)

This method is deprecated, please use monodomain_conductance() instead, but note the difference in the arguments:

calculate_conductance(r, sx, chi, dx)

is equivalent to:

monodomain_conductance(1 / chi, r, sx, dx, 1)
conductance()

Returns the cell-to-cell conductance(s) used in this simulation.

The returned value will be a single float for 1d simulations, a tuple (gx, gy) for 2d simulations, and None if conductances were set with set_conductance_field() or set_connections().

If diffusion is disabled, a call to this method will raise a RuntimeError.

default_state(x=None, y=None)

Returns the current default simulation state as a list of len(state) * n_total_cells floating point values, where n_total_cells is the total number of cells.

If the optional arguments x and y specify a valid cell index a single cell’s state is returned. For example state(4) can be used with a 1d simulation, while state(4, 2) is a valid index in the 2d case.

For 2d simulations, the list is indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

find_nan(log, watch_var=None, safe_range=None, return_log=False)

Searches for the origin of a bad value (NaN or inf) in a data log generated by this simulation.

Arguments:

log

A myokit.DataLog from this simulation. The log must contain the state of each cell and all bound variables. The bad value can occur at any point in time except the first.

watch_var

To aid in diagnosis, a variable can be selected as watch_var and a safe_range can be specified. With this option, the function will find and report either the first bad value or the first time the watched variable left the safe range, whatever came first. The watched variable must be a state variable.

safe_range

The safe range to check the watch_var against. The safe range should be specified as (lower, upper) where both bounds are assumed to be in the safe range.

return_log

If set to True, a log containing the final points before the error occurred will be returned, containing all variables for all cells.

Returns a tuple containing the following six values (or seven values if return_log=True):

time

The time the first bad value was found.

icell

The index of the cell in which the bad value was found, as a tuple e.g. (3, ) in a 1d simulation or (15, 12) in a 2d simulation.

variable

The name of the variable that was detected to be NaN or inf.

value

The bad value detected)

states

The state at time time and, if possible, up to 3 earlier states. Here, states[0] points to the current state, state[1] is the previous state and so on. Each state is represented as a list of values.

bounds

The bound variables corresponding to the returned states. For diffusion_current, this will contain the specific cell’s current, for time and pace the global values are shown. Each entry in bounds is a dictionary from variable names to values.

log

A myokit.DataLog with all variables for the same points as states and bounds. This will only be included if return_log=True.

is2d()

Deprecated alias of is_2d().

is_2d()

Returns True if and only if this is a 2d simulation.

is_paced(x, y=None)

Returns True if and only if the cell at index x (or index (x, y) in 2d simulations) will be paced during simulations.

If diffusion is disabled, a call to this method will raise a RuntimeError.

monodomain_conductance(chi, k, D, dx, A=1)

Calculates conductance values g based on monodomain parameters, assuming that membrane size, capacitance, and cell-to-cell conductance do not vary spatially.

In this simulation, the membrane potential V of a cell is assumed to be given by:

dV/dt = -1/Cm (I_ion + I_stim + I_diff)

where Cm is membrane capacitance and the diffusion current is defined by:

I_diff[i] = sum[g[i,j] * (V[i] - V[j])]

in which the sum is over all neighbors j of cell i.

Alternatively, with capacitance and currents normalized to membrane area, we can write:

dV/dt = -1/cm (i_ion + i_stim + i_diff)

In the monodomain model (see e.g. [4]), this diffusion current per unit area is defined as:

i_diff = -(1 / chi) (k / (k + 1)) ∇ * (D∇V)

(see the argument list below for the meaning of the variables).

This can be equated to Myokit’s diffusion current, but only if we assume zero-flux boundary conditions, a regularly spaced grid, and no spatial heterogeneity in D (or g).

With these assumptions, we can use finite differences to find:

g_bar = (1 / chi) * (k / (k + 1)) * D * (1 / dx^2)

where g_bar is the cell-to-cell conductance, but normalized with respect to unit membrane area. For models with currents normalized to area this is unproblematic, but to convert to models with unnormalized currents this means we have added the further assumption that each node contains some fixed amount of membrane, determined by an area A:

g = (1 / chi) * (k / (k + 1)) * D * (1 / dx^2) * A

This equation can also be applied in two dimensions, but only if we assume that the conductivity matrix is diagonal, in which case:

gx = (1 / chi) * (k / (k + 1)) * Dx * (1 / dx^2) * A
gy = (1 / chi) * (k / (k + 1)) * Dy * (1 / dy^2) * A

This method uses the above equation to calculate and return a conductance value from the parameters used in monodomain simulations.

Arguments:

chi

The surface area of the membrane per unit volume (in units area per volume).

k

The intra- to extracellular conductivity ratio.

D

The intracellular conductivity (in siemens per unit length).

dx

The distance between grid points.

A

The area of the cell membrane.

[4] Computer Model of Excitation and Recovery in the Anisotropic Myocardium I. Rectangular and Cubic Arrays of Excitable Elements. Leon & Horacek (1991) Journal of electrocardiology

neighbors(x, y=None)

Returns a list of indices specifying the neighbors of the cell at index x (or index (x, y) for 2d simulations).

Indices are given either as integers (1d or arbitrary geometry) or as tuples (2d).

If diffusion is disabled, a call to this method will raise a RuntimeError.

neighbours(x, y=None)

Alias of neighbors().

pre(duration, report_nan=True, progress=None, msg='Pre-pacing SimulationOpenCL')

This method can be used to perform an unlogged simulation, typically to pre-pace to a (semi-)stable orbit.

After running this method

  • The simulation time is not affected

  • The current state and the default state are updated to the final state reached in the simulation.

Calls to reset() after using pre() will revert the simulation to this new default state.

If numerical errors during the simulation lead to NaNs appearing in the result, the find_nan method will be used to pinpoint their location. Next, a call to the model’s rhs will be evaluated in python using checks for numerical errors to pinpoint the offending equation. The results of these operations will be written to stdout. To disable this feature, set report_nan=False.

To obtain feedback on the simulation progress, an object implementing the myokit.ProgressReporter interface can be passed in. passed in as progress. An optional description of the current simulation to use in the ProgressReporter can be passed in as msg.

remove_field(var)

Removes any field set for the given variable.

reset()

Resets the simulations:

  • The time variable is set to 0

  • The current state is set to the default state (either the model’s initial state or the last state reached using pre())

run(duration, log=None, log_interval=1.0, report_nan=True, progress=None, msg='Running SimulationOpenCL')

Runs a simulation and returns the logged results. Running a simulation has the following effects:

  • The internal state is updated to the last state in the simulation.

  • The simulation’s time variable is updated to reflect the time elapsed during the simulation.

The number of time units to simulate can be set with duration.

The variables to log can be indicated using the log argument. There are several options for its value:

  • None (default), to log all states

  • An integer flag or a combination of flags. Options: myokit.LOG_NONE, myokit.LOG_STATE, myokit.LOG_BOUND, myokit.LOG_INTER or myokit.LOG_ALL.

  • A list of qnames or variable objects

  • A myokit.DataLog object or another dictionary of

    qname : list mappings.

For more details on the log argument, see the function myokit.prepare_log().

Variables that vary from cell to cell will be logged with a prefix indicating the cell index. For example, when using:

s = SimulationOpenCL(m, p, ncells=256)
d = s.run(1000, log=['engine.time', 'membrane.V']

where engine.time is bound to time and membrane.V is the membrane potential variable, the resulting log will contain the following variables:

{
    'engine.time'  : [...],
    '0.membrane.V' : [...],
    '1.membrane.V' : [...],
    '2.membrane.V' : [...],
}

Alternatively, you can specify variables exactly:

d = s.run(1000, log=['engine.time', '0.membrane.V']

For 2d simulations, the naming scheme x.y.name is used, for example 0.0.membrane.V.

A log entry will be made every time at least log_interval time units have passed. No guarantee is given about the exact time log entries will be made, but the value of any logged time variable is guaranteed to be accurate.

If numerical errors during the simulation lead to NaNs appearing in the result, the find_nan method will be used to pinpoint their location. Next, a call to the model’s rhs will be evaluated in python using checks for numerical errors to pinpoint the offending equation. The results of these operations will be written to stdout. To disable this feature, set report_nan=False.

To obtain feedback on the simulation progress, an object implementing the myokit.ProgressReporter interface can be passed in. passed in as progress. An optional description of the current simulation to use in the ProgressReporter can be passed in as msg.

set_conductance(gx=10, gy=5)

Sets the cell-to-cell conductance used in this simulation.

For 1d simulations, only gx will be used and the argument gy can be omitted. For 2d simulations both arguments should be set.

The diffusion current is calculated as:

i = gx * ((V - V_xnext) - (V_xlast - V))
  + gy * ((V - V_ynext) - (V_ylast - V))

Where the second term gy * ... is only used for 2d simulations. At the boundaries, where either V_ilast or V_inext is unavailable, the value of V is substituted, causing the term to go to zero.

For a model with currents in [uA/uF] and voltage in [mV], gx` and gy have the unit [mS/uF].

Calling set_conductance will delete any conductances previously set with set_conductance_field() or set_connections().

If diffusion is disabled, a call to this method will raise a RuntimeError.

set_conductance_field(gx, gy=None)

Sets the cell-to-cell conductances used in this simulation, using lists of conductances.

For 1d simulations, the argument gx should be a sequence of nx - 1 non-negative floats (where nx is the number of cells), and gy should be None. For 2d simulations gx should be convertible to a numpy array with shape (ny, nx - 1) and gy should be convertible to a numpy array with shape (ny - 1, nx).

Calling set_conductance_field will delete any conductances previously set with set_conductance() or set_connections().

If diffusion is disabled, a call to this method will raise a RuntimeError.

set_connections(connections)

Adds a list of connections between cells, allowing the creation of arbitrary geometries.

The connections should be given as a list of tuples (cell_1, cell_2, conductance).

Connections are only supported in 1d mode – even though the resulting geometry may represent a shape in an arbitrary number of dimensions.

Calling set_connections will override any conductances previously set with set_conductance() or set_conductance_field().

If diffusion is disabled, a call to this method will raise a RuntimeError.

set_constant(var, value)

Changes a model constant. Only literal constants (constants not dependent on any other variable) can be changed.

The constant var can be given as a Variable or a string containing a variable qname. The value should be given as a float.

Note that any scalar fields set for the same variable will overwrite this value without warning.

set_default_state(state, x=None, y=None)

Changes this simulation’s default state.

This can be used in three different ways:

  1. When called with an argument state of size n_states and x=None the given state will be set as the new default state of all cells in the simulation.

  2. Called with an argument state of size n_states and x, y equal to a valid cell index, this method will update only the selected cell’s default state.

  3. Finally, when called with a state of size n_states * n_total_cells, the method will treat state as a concatenation of default state vectors for each cell. For 2d simulations, the list should be indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

set_field(var, values)

Can be used to replace a model constant with a scalar field.

The argument var must specify a variable from the simulation’s model. The field itself is given as values, which must have the dimensions (ny, nx). Multiple fields can be added, depending on the memory available on the device. If a field is added for a variable already associated with a field, the old data will be overwritten.

With diffusion currents enabled, this method can let you simulate heterogeneous tissue properties. With diffusion disabled, it can be used to investigate the effects of changing a parameter through the parallel simulation of several cells.

set_paced_cell_list(cells)

Selects the cells to be paced using a list of cell indices. In 1d simulations a cell index is an integer i, in 2d simulations cell indices are specified as tuples (i, j).

For large numbers of cells, this method becomes very inefficient. In these cases it may be better to use a rectangular pacing area set using set_paced_cells().

If diffusion is disabled, all cells are paced and a call to this method will raise a RuntimeError.

set_paced_cells(nx=5, ny=5, x=0, y=0)

Sets the number of cells that will receive a stimulus from the pacing protocol. For 1d simulations, the values ny and y will be ignored.

This method can only define rectangular pacing areas. To select an arbitrary set of cells, use set_paced_cell_list().

Arguments:

nx

The number of cells/nodes in the x-direction. If a negative number of cells is set the cells left of the offset (x) are stimulated.

ny

The number of cells/nodes in the y-direction. If a negative number of cells is set the cells left of the offset (x) are stimulated.

x

The offset of the pacing rectangle in the x-direction. If a negative offset is given the offset is calculated from right to left.

y

The offset of the pacing rectangle in the y-direction. If a negative offset is given the offset is calculated from bottom to top.

If diffusion is disabled, all cells are paced and a call to this method will raise a RuntimeError.

set_protocol(protocol=None)

Changes the pacing protocol used by this simulation.

set_state(state, x=None, y=None)

Changes the state of this simulation’s model.

This can be used in three different ways:

  1. When called with an argument state of size n_states and x=None the given state will be set as the new state of all cells in the simulation.

  2. Called with an argument state of size n_states and x, y equal to a valid cell index, this method will update only the selected cell’s state.

  3. Finally, when called with a state of size n_states * n_total_cells, the method will treat state as a concatenation of state vectors for each cell. For 2d simulations, the list should be indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

set_step_size(step_size=0.005)

Sets the step size used in the forward Euler solving routine.

set_time(time=0)

Sets the current simulation time.

shape()

Returns the shape of this simulation’s grid of cells as a tuple (ny, nx) for 2d simulations, or a single value nx for 1d simulations.

state(x=None, y=None)

Returns the current simulation state as a list of len(state) * n_total_cells floating point values, where n_total_cells is the total number of cells.

If the optional arguments x and y specify a valid cell index a single cell’s state is returned. For example state(4) can be used with a 1d simulation, while state(4, 2) is a valid index in the 2d case.

For 2d simulations, the list is indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

step_size()

Returns the current step size.

time()

Returns the current simulation time.

Fiber-Tissue Simulation

class myokit.FiberTissueSimulation(fiber_model, tissue_model, protocol=None, ncells_fiber=(128, 2), ncells_tissue=(128, 128), nx_paced=5, g_fiber=(9, 6), g_tissue=(9, 6), g_fiber_tissue=9, dt=0.005, precision=32, native_maths=False)

Runs a simulation of a fiber leading up to a rectangular piece of tissue.

Takes the following input arguments:

fiber_model

The model to simulate the fiber with.

tissue_model

The model to simulate the tissue with. Both models will be cloned when the simulation is created so that no changes to the given models will be made.

protocol

An optional pacing protocol, used to stimulate a number of cells at the start of the fiber.

ncells_fiber

The number of cells in the fiber (a tuple).

ncells_tissue

The number of cells in the tissue (a tuple).

nx_paced

The width (in cells) of the stimulus applied to the fiber. The fiber will be stimulated along its full height.

g_fiber

The cell to cell conductance in the fiber (a tuple).

g_tissue

The cell to cell conductance in the tissue (a tuple).

g_fiber_tissue

The fiber-cell to tissue-cell conductance at the junction (a scalar).

dt

The time step to use in the forward-Euler integration scheme.

precision

Can be set to myokit.SINGLE_PRECISION (default) or myokit.DOUBLE_PRECISION if the used device supports it.

native_maths

On some devices, selected functions (e.g. exp) can be made to run faster (but possibly less accurately) by setting native_maths=True.

The simulation provides the following inputs variables can bind to:

time (global)

The simulation time

pace (per-cell)

The pacing level, this is set if a protocol was passed in. Will be set to 0 if no protocol is provided.

diffusion_current (per-cell)

The current flowing from the cell to its neighbors. This will be positive when the cell is acting as a source, negative when it is acting as a sink. Will be set to 0 if no connections are made.

The variable time is set globally, meaning each cell uses the same value. The variables pace and diffusion_current have different values per cell.

The following labeled variables are required for this simulation to work:

membrane_potential

The variable representing the membrane potential.

The diffusion_current is calculated as:

i = gx * ((V - V_xnext) - (V_xlast - V))
  + gy * ((V - V_ynext) - (V_ylast - V))

At the boundaries, where either V_ilast or V_inext is unavailable, the value of V is substituted, causing the term to go to zero. The values of gx and gy can be set in the simulation’s constructor.

For a typical model with currents in [uA/uF] and voltage in [mV], gx` and gy have the unit [mS/uF].

Simulations maintain an internal state consisting of

  • the current simulation time

  • the current states

  • the default states

When a simulation is created, the simulation time is set to 0 and both the current and the default state are equal to the state of the given models, copied once for each cell. After each call to run() the time variable and current state are updated, so that each successive call to run continues where the previous simulation left off. A reset() method is provided that will set the time back to 0 and revert the current state to the default state. To change the time or state manually, use set_time() and set_fiber_state() and set_tissue_state().

A pre-pacing method pre() is provided that doesn’t affect the simulation time but will update the current and the default state. This allows you to pre-pace, run a simulation, reset to the pre-paced state, run another simulation etc.

The model passed to the simulation is cloned and stored internally, so changes to the original model object will not affect the simulation.

Models used with this simulation need to have independent components: it should be possible to evaluate the model’s equations one component at a time. A model’s suitability can be tested using :meth:` has_interdependent_components <myokit.Model.has_interdependent_components>`.

default_fiber_state(x=None, y=None)

Returns the current default simulation state for the fiber model as a list of len(state) * n_total_cells floating point values, where n_total_cells is the total number of fiber cells.

If the optional arguments x and y specify a valid cell index, a single cell’s state is returned.

The returned list is indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

default_tissue_state(x=None, y=None)

Returns the current default simulation state for the tissue model as a list of len(state) * n_total_cells floating point values, where n_total_cells is the total number of tissue cells.

If the optional arguments x and y specify a valid cell index, a single cell’s state is returned.

The returned list is indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

fiber_shape()

Returns the shape of this simulation’s grid of fiber cells as a tuple (ny, nx).

fiber_state(x=None, y=None)

Returns the current simulation state in the fiber as a list of len(state_fiber) * ncells_fiber floating point values.

If the optional arguments x and y specify a valid cell index a single fiber cell’s state is returned.

find_nan(logf, logt)

Searches for the origin of a NaN (or inf) in a set of simulation logs generated by this Simulation.

The logs must contain the state of each cell and all bound variables. The NaN can occur at any point in time except the first.

Returns a tuple (part, time, icell, variable, value, states, bound) where time is the time the first NaN was found and icell is the index of the cell in which it happened. The entry part is a string containing either “fiber” or “tissue”, indicating which part of the simulation triggered the error. The offending variable’s name is given as variable and its (illegal) value as value. The current state and, if available, any previous states are given in the list states. Here, states[0] points to the current state in the simulation part causing the error, state[1] is the previous state and so on. Similarly the values of the error causing model’s bound variables is given in bound.

pre(duration, report_nan=True, progress=None, msg='Pre-pacing FiberTissueSimulation')

This method can be used to perform an unlogged simulation, typically to pre-pace to a (semi-)stable orbit.

After running this method

  • The simulation time is not affected

  • The current state and the default state are updated to the final state reached in the simulation.

simulation to this new default state.

Calls to reset() after using pre() will revert the To obtain feedback on the simulation progress, an object implementing the myokit.ProgressReporter interface can be passed in. passed in as progress. An optional description of the current simulation to use in the ProgressReporter can be passed in as msg.

reset()

Resets the simulation:

  • The time variable is set to 0

  • The current state is set to the default state (either the model’s initial state or the last state reached using pre())

run(duration, logf=None, logt=None, log_interval=1.0, report_nan=True, progress=None, msg='Running FiberTissueSimulation')

Runs a simulation and returns the logged results as a tuple containing two myokit.DataLog objects.. Running a simulation has the following effects:

  • The internal state is updated to the last state in the simulation.

  • The simulation’s time variable is updated to reflect the time elapsed during the simulation.

The number of time units to simulate can be set with duration.

The variables to log can be indicated using the arguments logf and logt. There are several options for their values:

  • None (default), to log all states

  • An integer flag or a combination of flags. Options: myokit.LOG_NONE, myokit.LOG_STATE, myokit.LOG_BOUND, myokit.LOG_INTER, myokit.LOG_ALL.

  • A list of qnames or variable objects

  • A myokit.DataLog. In this case, new data will be appended to the existing log.

For more details on the log arguments, see the function myokit.prepare_log().

Any variables bound to “time” or “pace” will be logged globally, all others will be logged per cell. These variables will be prefixed with a single number indicating the cell index in the fiber, and with two numbers indicating the index in the tissue.

A log entry will be made every time at least log_interval time units have passed. No guarantee is given about the exact time log entries will be made, but the value of any logged time variable is guaranteed to be accurate.

Intermediary variables can be logged, but with one small drawback: for performance reasons the logged values of states and bound variables will always be one time step dt ahead of the intermediary variables. For example if running the simulation with a step size dt=0.001 the entry for a current IKr stored at t=1 will be IKr(0.999), while the entry for state V will be V(1). If exact intermediary variables are needed it’s best to log only states and bound variables and re-calculate the intermediary variables from these manually.

To obtain feedback on the simulation progress, an object implementing the myokit.ProgressReporter interface can be passed in. passed in as progress. An optional description of the current simulation to use in the ProgressReporter can be passed in as msg.

set_default_fiber_state(state, x=None, y=None)

Changes this simulation’s default state for the fiber model.

This can be used in three different ways:

  1. When called with an argument state of size n_states and x=None the given state will be set as the new default state of all fiber cells in the simulation.

  2. Called with an argument state of size n_states and x, y equal to a valid cell index, this method will update only the selected fiber cell’s default state.

  3. Finally, when called with a state of size n_states * n_total_cells, the method will treat state as a concatenation of default state vectors for each fiber cell. For 2d simulations, the list should be indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

set_default_tissue_state(state, x=None, y=None)

Changes this simulation’s default state for the tissue model.

This can be used in three different ways:

  1. When called with an argument state of size n_states and x=None the given state will be set as the new default state of all tissue cells in the simulation.

  2. Called with an argument state of size n_states and x, y equal to a valid cell index, this method will update only the selected tissue cell’s default state.

  3. Finally, when called with a state of size n_states * n_total_cells, the method will treat state as a concatenation of default state vectors for each tissue cell. For 2d simulations, the list should be indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

set_fiber_state(state, x=None, y=None)

Changes the state of this simulation’s fiber model.

This can be used in three different ways:

  1. When called with an argument state of size n_states and x=None the given state will be set as the new state of all fiber cells in the simulation.

  2. Called with an argument state of size n_states and x, y equal to a valid cell index, this method will update only the selected fiber cell’s state.

  3. Finally, when called with a state of size n_states * n_total_cells, the method will treat state as a concatenation of state vectors for each fiber cell. For 2d simulations, the list should be indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

set_protocol(protocol=None)

Changes the pacing protocol used by this simulation.

set_step_size(step_size=0.005)

Sets the solver step size.

set_time(time=0)

Sets the current simulation time.

set_tissue_state(state, x=None, y=None)

Changes the state of this simulation’s tissue model.

This can be used in three different ways:

  1. When called with an argument state of size n_states and x=None the given state will be set as the new state of all tissue cells in the simulation.

  2. Called with an argument state of size n_states and x, y equal to a valid cell index, this method will update only the selected tissue cell’s state.

  3. Finally, when called with a state of size n_states * n_total_cells, the method will treat state as a concatenation of state vectors for each tissue cell. For 2d simulations, the list should be indexed so that x changes first (the first nx entries have y = 0, the second nx entries have y = 1 and so on).

step_size(step_size=0.005)

Returns the solver step size.

time()

Returns the current simulation time

tissue_shape()

Returns the shape of this simulation’s grid of tissue cells as a tuple (ny, nx).

tissue_state(x=None, y=None)

Returns the current simulation state in the tissue as a list of len(state_tissue) * ncells_tissue floating point values.

If the optional arguments x and y specify a valid cell index a single cell’s state is returned.

OpenCL utility classes

class myokit.OpenCL

Tests for OpenCL support and can return information about opencl simulations.

static available()

Returns True if OpenCL support has been detected on this sytem and at least one platform and device were detected.

static current_info(formatted=False)

Returns a myokit.OpenCLPlatformInfo object for the platform and device selected by the user, or chosen as default.

If formatted=True is set, a formatted version of the information is returned instead.

static info(formatted=False)

Queries the OpenCL installation for the available platforms and devices and returns a myokit.OpenCLInfo object.

If formatted=True is set, a formatted version of the information is returned instead.

static load_selection()

Loads a platform/device selection from disk and returns a tuple (platform, device). Each entry in the tuple is either a string with the platform/device name, or None if no preference was set.

static load_selection_bytes()

Loads a platform/device selection from disk and returns a tuple (platform, device). Each entry in the tuple is either a string with the platform/device name, or None if no preference was set.

static save_selection(platform=None, device=None)

Stores a platform/device selection to disk.

Both platform and device are identified by their names.

static selection_info()

Returns a list of platform/device combinations along with information allowing the user to select one.

The returned output is a list of tuples, where each tuple has the form (platform_name, device_name, specs).

A preferred device can be selected by passing one of the returned platform_name, device_name combinations to OpenCL.set_preferred_device().

static supported()

Returns True if OpenCL support has been detected on this system.

class myokit.OpenCLInfo(mcl_info=[])

Represents information about the available OpenCL platforms and devices.

Each OpenCLInfo object has a property platforms, containing a tuple of OpenCLPlatformInfo objects.

OpenCLInfo objects are created and returned by myokit.OpenCL.

format()

Returns a formatted string version of this object’s information.

class myokit.OpenCLPlatformInfo(platform)

Represents information about an OpenCL platform.

An OpenCLPlatformInfo object has the following properties:

name (string)

This platform’s name.

vendor (string)

The vendor of this platform.

version (string)

The OpenCL version supported by this platform.

profile (string)

The supported OpenCL profile of this platform.

extensions (string)

The available OpenCL extensions on this platform.

devices (tuple)

A tuple of device information objects for the devices available on this platform. This field may be None, in which case device will be set instead.

device (OpenCLDeviceInfo)

An information objects for the device selected by the user, or chosen as the default device. This field may be None, in which case devices will be set instead.

OpenCLPlatformInfo objects are created as part of a OpenCLInfo objects, as returned by most OpenCL enabled parts of Myokit.

format()

Returns a formatted string version of this object’s information.

class myokit.OpenCLDeviceInfo(device)

Represents information about an OpenCL device.

An OpenCLDeviceInfo object has the following properties:

name (string)

This device’s name.

vendor (string)

This device’s vendor.

version (string)

The OpenCL version supported by this device.

driver (string)

The driver version for this device.

clock (int)

This device’s clock speed (in MHz).

globl (int)

The available global memory on this device (in bytes).

local (int)

The available local memory on this device (in bytes).

const (int)

The available constant memory on this device (in bytes).

units (int)

The number of computing units on this device.

param (int)

The maximum total size (in bytes) of arguments passed to the kernel. This limits the number of arguments a kernel can get.

groups (int)

The maximum work group size.

dimensions (int)

The maximum work item dimension.

items (tuple)

A tuple of ints specifying the maximum work item size in each dimension.

OpenCLDeviceInfo objects are created as part of a OpenCLInfo objects, as returned by most OpenCL enabled parts of Myokit.