Cable simulation¶
- class myokit.Simulation1d(model, protocol=None, ncells=50, rl=False)¶
Can run 1d cable simulations based on a
model
.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 at the start of the cable.
ncells
The number of cells in the cable
rl
Use Rush-Larsen updates instead of forward-Euler for any Hodgkin-Huxley gating variables (default=False).
This 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
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 variablespace
anddiffusion_current
have different values per cell. The number of paced cell can be set withset_paced_cells()
.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. Areset()
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, useset_time()
andset_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
diffusion_current
is calculated as:i = sum[g * (V - V_j)]
Where the sum is taken over all neighboring cells j (see [1]).
The resulting ODE system is solved using a forward Euler (FE) method with fixed step sizes. Smaller step sizes lead to more accurate results, and it is recommended any important results are double-checked by re-running with a reduced step size. Any states written in a Hodgkin-Huxley form can be updated using Rush-Larsen steps (see [2]), by setting
rl=True
. This often increases stability (allowing for larger step sizes) but can reduce 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
- conductance()¶
Returns the current conductance.
- default_state(icell=None)¶
Returns the default simulation state as a list of
len(state) * ncells
floating point values. If the optional argumenticell
is set to a valid cell index only the state of that cell is returned.
- paced_cells()¶
Returns the number of cells that will receive a stimulus from the pacing protocol.
- pre(duration, progress=None, msg='Pre-pacing Simulation1d')¶
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 usingpre()
will revert the simulation to this new default state.To obtain feedback on the simulation progress, an object implementing the
myokit.ProgressReporter
interface can be passed in. passed in asprogress
. 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 state is set to the default state
- run(duration, log=None, log_interval=1.0, progress=None, msg='Running Simulation1d')¶
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 statesAn integer flag or a combination of flags. Options:
myokit.LOG_NONE
,myokit.LOG_STATE
,myokit.LOG_INTER
,myokit.LOG_BOUND
,myokit.LOG_ALL
.A list of qnames or variable objects
A
myokit.DataLog
object.
For more details on the
log
argument, see the functionmyokit.prepare_log()
.Any variables bound to “time” will be logged globally, all others will be logged per cell. These variables will be prefixed with a number indicating the cell index. For example, when using:
s = Simulation1d(m, p, ncells=3) 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' : [...], }
A log entry is created every time at least
log_interval
time units have passed.To obtain feedback on the simulation progress, an object implementing the
myokit.ProgressReporter
interface can be passed in. passed in asprogress
. An optional description of the current simulation to use in the ProgressReporter can be passed in as msg.
- set_conductance(g=10)¶
Changes the cell-to-cell conductance.
- set_default_state(state, icell=None)¶
Changes this simulation’s default state.
This can be used in three different ways:
When called with an argument
state
of sizen_states
andi_cell=None
the given state will be set as the new state of all cells in the simulation.Called with an argument
state
of size n_states andi_cell
equal to a valid cell index, this method will update only the selected cell’s state.Finally, when called with a
state
of sizen_states * ncells
the method will treatstate
as a concatenation of state vectors for each cell.
- set_paced_cells(n=5)¶
Sets the number of cells that will receive a stimulus from the pacing protocol.
- set_protocol(protocol=None)¶
Changes the pacing protocol used by this simulation.
- set_state(state, icell=None)¶
Changes the state of this simulation’s model.
This can be used in three different ways:
When called with an argument
state
of sizenstates
andicell=None
the given state will be set as the new state of all cells in the simulation.Called with an argument
state
of sizenstates
andicell
equal to a valid cell index, this method will update only the selected cell’s state.Finally, when called with a
state
of sizenstates * ncells
the method will treatstate
as a concatenation of state vectors for each cell.
- set_step_size(step_size=0.005)¶
Sets the solver step size. In some cases, the solver will take a slightly smaller step size, either to arrive exactly at the start/end of a pacing event or to arrive exactly at the end of a simulation.
- set_time(time=0)¶
Sets the current simulation time.
- state(icell=None)¶
Returns the current simulation state as a list of
len(state) * ncells
floating point values. If the optional argumenticell
is set to a valid cell index only the state of that cell is returned.
- time()¶
Returns the current simulation time.