FEniCSx-Beat#

class beat.ECGRecovery(v: dolfinx.fem.function.Function, sigma_b: float | dolfinx.fem.function.Constant = 1.0, C_m: float | dolfinx.fem.function.Constant = 1.0, dx: ufl.measure.Measure | None = None, M: float = 1.0, petsc_options: dict[str, typing.Any] = <factory>)[source]#
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, dx: Measure | None = None, **kwargs)[source]#

Solve

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

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

Parameters:

dt (float) – Time step size

Returns:

The variational form and the precondition

Return type:

tuple[ufl.Form, ufl.Form]

class beat.MonodomainSplittingSolver(pde: beat.monodomain_model.MonodomainModel, ode: beat.monodomain_solver.ODESolver, theta: float = 1.0)[source]#
class beat.Stimulus(expr, dZ, marker)[source]#
dZ: Measure#

Alias for field number 1

expr: Expr#

Alias for field number 0

marker: int | None#

Alias for field number 2

class beat.base_model.BaseModel(time: Constant, mesh: Mesh, dx: Measure | None = None, params: dict[str, Any] | None = None, I_s: Stimulus | Sequence[Stimulus] | Coefficient | None = None, jit_options: dict[str, Any] | None = None, form_compiler_options: dict[str, Any] | None = None, petsc_options: dict[str, Any] | None = None)[source]#

Base class for models.

Parameters:
  • time (dolfinx.fem.Constant) – The current time

  • mesh (dolfinx.mesh.Mesh) – The mesh

  • dx (ufl.Measure, optional) – The measure for the spatial domain, by default None

  • params (dict, optional) – Parameters for the model, by default None

  • I_s (Stimulus | Sequence[Stimulus] | ufl.Coefficient, optional) – The stimulus, by default None

  • jit_options (dict, optional) – JIT options, by default None

  • form_compiler_options (dict, optional) – Form compiler options, by default None

  • petsc_options (dict, optional) – PETSc options, by default None

solve(interval: tuple[float, float], dt: float | None = None) Results[source]#

Solve on the given time interval.

Parameters:
  • interval (tuple[float, float]) – The time interval (T0, T)

  • dt (float, optional) – The time step, by default None

Returns:

The results of the solution

Return type:

Results

step(interval)[source]#

Perform a single time step.

Parameters:

interval (tuple[float, float]) – The time interval (T0, T)

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

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

Parameters:

dt (Expr | float) – The time step

Returns:

The variational form and the precondition

Return type:

tuple[ufl.Form, ufl.Form]

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, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)[source]#
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, dx: Measure | None = None, **kwargs)[source]#

Solve

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

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

Parameters:

dt (float) – Time step size

Returns:

The variational form and the precondition

Return type:

tuple[ufl.Form, ufl.Form]

class beat.monodomain_solver.MonodomainSplittingSolver(pde: beat.monodomain_model.MonodomainModel, ode: beat.monodomain_solver.ODESolver, theta: float = 1.0)[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: 'dolfinx.fem.Function', v_pde: 'dolfinx.fem.Function', markers: 'dolfinx.fem.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 dolfinx function to numpy array

to_dolfin() None[source]#

Assign values from numpy array to dolfinx function

class beat.odesolver.DolfinODESolver(v_ode: 'dolfinx.fem.Function', v_pde: 'dolfinx.fem.Function', init_states: 'npt.NDArray', parameters: 'npt.NDArray', fun: 'Callable', num_states: 'int', v_index: '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')[source]#
beat.utils.expand_layer(V: FunctionSpace, ft: MeshTags, endo_marker: int, epi_marker: int, endo_size: float, epi_size: float, output_mid_marker: int = 0, output_endo_marker: int = 1, output_epi_marker: int = 2) Function[source]#

Expand the endo and epi markers to the rest of the mesh with a given size

Parameters:
  • V (dolfinx.fem.FunctionSpace) – The function space for your return function

  • ft (dolfinx.mesh.MeshTags) – The facet tags

  • endo_marker (int) – The endocardium marker

  • epi_marker (int) – The epicardium marker

  • endo_size (float) – The endocardium size

  • epi_size (float) – The epicardium size

  • output_mid_marker (int, optional) – The marker to be set for the midwall, by default 0

  • output_endo_marker (int, optional) – The marker to be set for the endocardium, by default 1

  • output_epi_marker (int, optional) – The marker to be set for the epicardium, by default 2

Returns:

The expanded markers as a function in V

Return type:

dolfinx.fem.Function

beat.utils.expand_layer_biv(V: FunctionSpace, ft: MeshTags, endo_lv_marker: int, endo_rv_marker: int, epi_marker: int, endo_size: float, epi_size: float, output_mid_marker: int = 0, output_endo_marker: int = 1, output_epi_marker: int = 2) Function[source]#

Expand the endo and epi markers to the rest of the mesh with a given size

Parameters:
  • V (dolfinx.fem.FunctionSpace) – The function space for your return function

  • ft (dolfinx.mesh.MeshTags) – The facet tags

  • endo_lv_marker (int) – The LV endocardium marker

  • endo_rv_marker (int) – The LV endocardium marker

  • epi_marker (int) – The epicardium marker

  • endo_size (float) – The endocardium size

  • epi_size (float) – The epicardium size

  • output_mid_marker (int, optional) – The marker to be set for the midwall, by default 0

  • output_endo_marker (int, optional) – The marker to be set for the endocardium, by default 1

  • output_epi_marker (int, optional) – The marker to be set for the epicardium, by default 2

Returns:

The expanded markers as a function in V

Return type:

dolfinx.fem.Function

beat.utils.local_project(v: Function, V: FunctionSpace, u: Function | None = None) Function | None[source]#

Element-wise projection using LocalSolver

Parameters:
  • v (dolfinx.fem.Function) – Function to be projected

  • V (dolfinx.fem.FunctionSpace) – Function space to project into

  • u (dolfinx.fem.Function | None, optional) – Optional function to save the projected function, by default None

Returns:

The projected function

Return type:

dolfinx.fem.Function | None

beat.utils.parse_element(space_string: str, mesh: Mesh, dim: int) _ElementBase[source]#

Parse a string representation of a basix element family

beat.utils.space_from_string(space_string: str, mesh: Mesh, dim: int = 1) functionspace[source]#

Constructed a finite elements space from a string representation of the space

Parameters:
  • space_string (str) – A string on the form {family}_{degree} which determines the space. Example ‘Lagrange_1’.

  • mesh (df.Mesh) – The mesh

  • dim (int) – 1 for scalar space, 3 for vector space.

Returns:

The function space

Return type:

df.FunctionSpace