API Reference#

Top-level package for cbcbeatx.

class cbcbeatx.Markerwise(objects: list[Expr], keys: list[int], marker: MeshTags)[source]#

A container class representing an object defined by a number of objects combined with a mesh tag defining mesh domains and a map from each object to the subdomain.

Parameters:
  • objects (list[ufl.core.expr.E[xr]) – The forcing terms as ufl expressions

  • keys (list[int]) – The cell-tag integer associated with each forcing term

  • marker (dolfinx.mesh.MeshTags) – The cell tags

Examples

Given (g0, g1), (2, 5) and cell_markers, let

\[\begin{split}g =\begin{cases} g0 \text{ on domains marked by 2 in } cell\_markers \\ g1 \text{ on domains marked by 5 in } cell\_markers \end{cases}\end{split}\]
g = Markerwise((g0, g1), (2, 5), markers)
property marker: MeshTags#

The cell marker

class cbcbeatx.MonodomainSolver(mesh: Mesh, M_i: Expr, I_s: tuple[Expr, Markerwise] | None = None, v0: Expr | None = None, time: Constant | None = None, params: dict | None = None)[source]#

Solve the (pure) monodomain equations on the form: Find the transmembrane potential \(v = v(x, t)\) such that

\[\frac{\partial v}{\partial t} - \nabla \cdot ( M_i \nabla v) = I_s,\]

where \(M_i\) is the intracellular cardiac conductivity tensor; \(I_s\) is prescribed input. In addition, initial conditions are given for \(v\):

\[v(x, 0) = v_0\]

This solver assumes pure homogeneous Neumann boundary conditions for \(v\).

This solver is based on a \(\theta\)-scheme discretization in time and N-th order Lagrange space elements in space.

Parameters:
  • mesh – The spatial domain (mesh).

  • M_i – The intracellular conductivity tensor (as an UFL expression).

  • time – A constant holding the current time. If None is given, time is created for you, initialized to zero. This can be used as an internal variable in the RHS term \(I_s\) to get a time-dependent input.

  • I_s – A typically time-dependent external stimulus given as a dict, with domain markers as the key and a dolfinx.fem.Expression as values. NB: it is assumed that the time dependence of I_s is encoded via the ‘time’ Constant.

  • v0 – Initial condition for v.

  • params

    Parameter dictionary specifying problem specific parameters.

    "polynomial_degree": The degree of the transmembral potential.

    "theta": Degree used in theta scheme.

    "jit_options": Dictionary with options for Just in time compilation. See dolfinx/jit.py for more information.

    "form_compiler_options": Dictionary with options for the FFCx form compiler. Call python3 -m ffcx --help for all available options.

    "petsc_options": Dictionary containing options for the PETSc KSP solver. See the PETSc documentation for more information.

    "use_custom_preconditioner": True/False Use \(\int_\Omega v\cdot w + \frac{\Delta t}{2} (M_i \nabla v)\cdot \nabla w~\mathrm{d}x\) as preconditioner

    dt: Initial time step

Examples

import dolfinx
import ufl
from mpi4py import MPI
from petsc4py import PETSc
from cbcbeatx import MonodomainSolver
N = 10
mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, N, N)
M_i = ufl.as_tensor(((0.2, 0), (0., 0.3)))
time = dolfinx.fem.Constant(mesh, PETSc.ScalarType(0.))
x = ufl.SpatialCoordinate(mesh)
I_s = x[0] * ufl.cos(x[1]) * time
v_0 = ufl.sin(2*ufl.pi*x[0])
petsc_options = {"ksp_type": "preonly", "pc_type": "lu"}
params = {"polynomial_degree": 2, "theta": 0.5, "petsc_options": petsc_options,
"use_custom_preconditioner": True}
solver = MonodomainSolver(mesh, M_i, time, I_s, v_0, params)
static default_parameters() dict[source]#

Get the default parameters for the class

Returns:

The default parameters

Return type:

dict

solution_fields() tuple[Function, Function][source]#

Return a tuple of Functions \((u_{n-1}, u_n)\) where \(u_{n-1}\) is the solution from the previous time step. \(u_n\) the solution at the current time step

Returns:

The functions

Return type:

tuple[dolfinx.fem.Function, dolfinx.fem.Function]

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

Solve the problem on a time given interval \([T_0, T_1]\)

Parameters:
  • interval (tuple[float, float]) – The time interval \([T_0, T_1]\)

  • dt (float, optional) – The time step \(\Delta t\). Defaults to None, which corresponds to \(\Delta t = T_1-T_0\).

Yields:

step_output

An iterator solving for each time step. Each element of the iterator is

a tuple describing the time step, and a tuple of the previous and current solution.

step(interval: tuple[float, float])[source]#

Solve the problem on a time given interval \([T_0, T_1]\)

Parameters:
  • interval (tuple[float, float]) – The time interval \([T_0, T_1]\)

  • dt (float, optional) – The time step \(\Delta t\). Defaults to None, which corresponds to \(\Delta t = T_1-T_0\).

property time: float64#

Return current time used in solver

cbcbeatx.rhs_with_markerwise_field(V: FunctionSpace, g: Expr | Markerwise | None) tuple[Measure, Form][source]#

Create the ufl-form \(G=\int_\Omega g v~\mathrm{d}x\) where g can be:

  1. A set of ufl expressions over subdomains (Markerwise)

  2. A single ufl expression over the whole domain (ufl.core.expr.Expr)

  3. If no expression is supplied, return a zero integral (will be optimized away later)

Parameters:
  • V – The function space to extract the test function v from

  • g – The forcing term

Returns:

A tuple of the integration measure dx over the domain (without subdomain id set) and the source form G