FEniCS-Beat#

class beat.BidomainModel(mesh, time, M_i, M_e, I_s: Stimulus | Coefficient | None = None, I_a=None, v_=None, params=None)[source]#
static default_parameters()[source]#

Initialize and return a set of default parameters

Returns

A set of parameters (dolfin.Parameters)

To inspect all the default parameters, do:

info(BidomainSolver.default_parameters(), True)
variational_forms(k_n: Expr | float) tuple[Form, Form][source]#

Create the variational forms corresponding to the given discretization of the given system of equations.

Arguments
k_n (ufl.Expr or float)

The time step

Returns

(lhs, rhs) (tuple of ufl.Form)

class beat.Geometry(mesh, ffun, markers, f0, s0, n0)[source]#
f0: Constant | Function | None#

Alias for field number 3

ffun: MeshFunction | None#

Alias for field number 1

markers: dict[str, tuple[int, int]] | None#

Alias for field number 2

mesh: Mesh#

Alias for field number 0

n0: Constant | Function | None#

Alias for field number 5

s0: Constant | Function | None#

Alias for field number 4

class beat.MonodomainModel(time: Constant, mesh: Mesh, M: Coefficient | float, I_s: Stimulus | Sequence[Stimulus] | Coefficient | None = None, params=None, C_m: float = 1.0)[source]#

Solve

\[\frac{\partial V}{\partial t} - \nabla \cdot \left( M \nabla V \right) - I_{\mathrm{stim}} = 0\]
variational_forms(k_n: Expr | float) tuple[Form, Form][source]#

Create the variational forms corresponding to the given discretization of the given system of equations.

Arguments
k_n (ufl.Expr or float)

The time step

Returns

(lhs, rhs, prec) (tuple of ufl.Form)

class beat.MonodomainSplittingSolver(pde: beat.monodomain_model.MonodomainModel, ode: beat.monodomain_solver.ODESolver, theta: float = 0.5)[source]#
class beat.base_model.BaseModel(time: Constant, mesh: Mesh, params: dict[str, Any] | None = None, I_s: Stimulus | Sequence[Stimulus] | Coefficient | None = None)[source]#
solve(interval: tuple[float, float], dt: float | None = None)[source]#

Solve the discretization on a given time interval (t0, t1) with a given timestep dt and return generator for a tuple of the interval and the current solution.

Arguments
interval (tuple)

The time interval for the solve given by (t0, t1)

dt (int, optional)

The timestep for the solve. Defaults to length of interval

Returns

(timestep, solution_field) via (genexpr)

Example of usage:

# Create generator
solutions = solver.solve((0.0, 1.0), 0.1)

# Iterate over generator (computes solutions as you go)
for (interval, solution_fields) in solutions:
  (t0, t1) = interval
  v_, v = solution_fields
  # do something with the solutions
step(interval)[source]#

Solve on the given time step (t0, t1).

abstract variational_forms(k_n: Expr | float) tuple[Form, Form][source]#

Create the variational forms corresponding to the given discretization of the given system of equations.

class beat.base_model.Results(state, status)[source]#
state: Function#

Alias for field number 0

status: Status#

Alias for field number 1

class beat.base_model.Status(value)[source]#

An enumeration.

class beat.bidomain_model.BidomainModel(mesh, time, M_i, M_e, I_s: Stimulus | Coefficient | None = None, I_a=None, v_=None, params=None)[source]#
static default_parameters()[source]#

Initialize and return a set of default parameters

Returns

A set of parameters (dolfin.Parameters)

To inspect all the default parameters, do:

info(BidomainSolver.default_parameters(), True)
variational_forms(k_n: Expr | float) tuple[Form, Form][source]#

Create the variational forms corresponding to the given discretization of the given system of equations.

Arguments
k_n (ufl.Expr or float)

The time step

Returns

(lhs, rhs) (tuple of ufl.Form)

class beat.ecg.Leads12(RA, LA, LL, RL, V1, V2, V3, V4, V5, V6)[source]#
property I: ndarray#

Voltage between the (positive) left arm (LA) electrode and right arm (RA) electrode

property II: ndarray#

Voltage between the (positive) left leg (LL) electrode and the right arm (RA) electrode

property III: ndarray#

Voltage between the (positive) left leg (LL) electrode and the left arm (LA) electrode

LA: ndarray#

Alias for field number 1

LL: ndarray#

Alias for field number 2

RA: ndarray#

Alias for field number 0

RL: ndarray | None#

Alias for field number 3

V1: ndarray | None#

Alias for field number 4

V2: ndarray | None#

Alias for field number 5

V3: ndarray | None#

Alias for field number 6

V4: ndarray | None#

Alias for field number 7

V5: ndarray | None#

Alias for field number 8

V6: ndarray | None#

Alias for field number 9

property Vw: ndarray#

Wilson’s central terminal

property aVF: ndarray#

Lead augmented vector foot (aVF) has the positive electrode on the left leg. The negative pole is a combination of the right arm electrode and the left arm electrode

property aVL: ndarray#

Lead augmented vector left (aVL) has the positive electrode on the left arm. The negative pole is a combination of the right arm electrode and the left leg electrode

property aVR: ndarray#

Lead augmented vector right (aVR) has the positive electrode on the right arm. The negative pole is a combination of the left arm electrode and the left leg electrode

class beat.monodomain_model.MonodomainModel(time: Constant, mesh: Mesh, M: Coefficient | float, I_s: Stimulus | Sequence[Stimulus] | Coefficient | None = None, params=None, C_m: float = 1.0)[source]#

Solve

\[\frac{\partial V}{\partial t} - \nabla \cdot \left( M \nabla V \right) - I_{\mathrm{stim}} = 0\]
variational_forms(k_n: Expr | float) tuple[Form, Form][source]#

Create the variational forms corresponding to the given discretization of the given system of equations.

Arguments
k_n (ufl.Expr or float)

The time step

Returns

(lhs, rhs, prec) (tuple of ufl.Form)

class beat.monodomain_solver.MonodomainSplittingSolver(pde: beat.monodomain_model.MonodomainModel, ode: beat.monodomain_solver.ODESolver, theta: float = 0.5)[source]#
class beat.monodomain_solver.ODESolver(*args, **kwargs)[source]#
class beat.odesolver.BaseDolfinODESolver[source]#
ode_to_pde() None[source]#

Projects v_ode (DG0, quadrature space, …) into v_pde (CG1)

pde_to_ode() None[source]#

Projects v_pde (CG1) into v_ode (DG0, quadrature space, …)

class beat.odesolver.DolfinMultiODESolver(v_ode: 'dolfin.Function', v_pde: 'dolfin.Function', markers: 'dolfin.Function', init_states: 'dict[int, npt.NDArray]', parameters: 'dict[int, npt.NDArray]', fun: 'dict[int, Callable]', num_states: 'dict[int, int]', v_index: 'dict[int, int]')[source]#
from_dolfin() None[source]#

Assign values from dolifn function to numpy array

to_dolfin() None[source]#

Assign values from numpy array to dolfin function

class beat.odesolver.DolfinODESolver(v_ode: 'dolfin.Function', v_pde: 'dolfin.Function', init_states: 'npt.NDArray', parameters: 'npt.NDArray', fun: 'Callable', num_states: 'int', v_index: 'int' = 0, missing_variables: 'npt.NDArray | None' = None, num_missing_variables: 'int' = 0)[source]#
from_dolfin() None[source]#

Assign values from dolfin function to numpy array

to_dolfin() None[source]#

Assign values from numpy array to dolfin function

class beat.odesolver.ODEResults(y, t)[source]#
t: ndarray[Any, dtype[float64]]#

Alias for field number 1

y: ndarray[Any, dtype[float64]]#

Alias for field number 0

class beat.odesolver.ODESystemSolver(fun: 'Callable', states: 'npt.NDArray', parameters: 'npt.NDArray', missing_variables: 'npt.NDArray | None' = None, _kwargs: 'dict[str, npt.NDArray]' = <factory>)[source]#