FEniCSx-LDRB#

ldrb.dolfinx_ldrb(mesh: Mesh, fiber_space: str = 'P_1', ffun: MeshTags | None = None, markers: dict[str, list[int]] | dict[str, tuple[int, ...]] | dict[str, int] | None = None, alpha_endo_lv: float = 40, alpha_epi_lv: float = -50, alpha_endo_rv: float | None = None, alpha_epi_rv: float | None = None, alpha_endo_sept: float | None = None, alpha_epi_sept: float | None = None, beta_endo_lv: float = -65, beta_epi_lv: float = 25, beta_endo_rv: float | None = None, beta_epi_rv: float | None = None, beta_endo_sept: float | None = None, beta_epi_sept: float | None = None, epi_only: bool = False, **kwargs) LDRBOutput[source]#

Create fiber, cross fibers and sheet directions

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

  • fiber_space (str) – A string on the form {family}_{degree} which determines for what space the fibers should be calculated for. If not provided, then a first order Lagrange space will be used, i.e Lagrange_1.

  • ffun (dolfinx.mesh.MeshTags) – A facet function containing markers for the boundaries.

  • markers (dict[str, int | list[int]] (optional)) – A dictionary with the markers for the different boundaries defined in the facet function or within the mesh itself. The following markers must be provided: ‘base’, ‘lv’, ‘epi, ‘rv’ (optional). If the markers are not provided the following default vales will be used: base = 10, rv = 20, lv = 30, epi = 40

  • save_markers (bool) – If true save markings of the geometry. This is nice if you want to see that the LV, RV and Septum are marked correctly.

  • angles (kwargs) –

    Keyword arguments with the fiber and sheet angles. It is possible to set different angles on the LV, RV and septum, however it either the RV or septum angles are not provided, then the angles on the LV will be used. The default values are taken from the original paper, namely

    \[\begin{split}\alpha_{\text{endo}} &= 40 \\ \alpha_{\text{epi}} &= -50 \\ \beta_{\text{endo}} &= -65 \\ \beta_{\text{epi}} &= 25\end{split}\]

    The following keyword arguments are possible:

    alpha_endo_lvscalar

    Fiber angle at the LV endocardium.

    alpha_epi_lvscalar

    Fiber angle at the LV epicardium.

    beta_endo_lvscalar

    Sheet angle at the LV endocardium.

    beta_epi_lvscalar

    Sheet angle at the LV epicardium.

    alpha_endo_rvscalar

    Fiber angle at the RV endocardium.

    alpha_epi_rvscalar

    Fiber angle at the RV epicardium.

    beta_endo_rvscalar

    Sheet angle at the RV endocardium.

    beta_epi_rvscalar

    Sheet angle at the RV epicardium.

    alpha_endo_septscalar

    Fiber angle at the septum endocardium.

    alpha_epi_septscalar

    Fiber angle at the septum epicardium.

    beta_endo_septscalar

    Sheet angle at the septum endocardium.

    beta_epi_septscalar

    Sheet angle at the septum epicardium.

  • epi_only (bool) – If true, then only the angles on the epicardial surface will be used. This is a hack to make it possible to compute purely longitudinal fibers

ldrb.project_gradients(mesh: Mesh, scalar_solutions: dict[str, Function], fiber_space: str = 'P_1') tuple[dict[str, ndarray], dict[str, Function]][source]#

Calculate the gradients using projections

Parameters:
  • mesh (dolfin.Mesh) – The mesh

  • fiber_space (str) – A string on the form {family}_{degree} which determines for what space the fibers should be calculated for.

  • scalar_solutions (dict) – A dictionary with the scalar solutions that you want to compute the gradients of.

ldrb.scalar_laplacians(mesh: Mesh, markers: dict[str, list[int]], ffun: MeshTags) dict[str, Function][source]#

Calculate the laplacians

Parameters:
  • mesh (dolfin.Mesh) – A dolfin mesh

  • markers (dict (optional)) – A dictionary with the markers for the different boundaries defined in the facet function or within the mesh itself. The following markers must be provided: ‘base’, ‘lv’, ‘epi, ‘rv’ (optional). If the markers are not provided the following default vales will be used: base = 10, rv = 20, lv = 30, epi = 40.

ldrb.space_from_string(space_string: str, mesh: Mesh, dim: int) 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.

ldrb.calculus.axis(u: ndarray, v: ndarray) ndarray[source]#

Construct the fiber orientation coordinate system.

Given two vectors \(u\) and \(v\) in the apicobasal and transmural direction respectively return a matrix that represents an orthonormal basis in the circumferential (first), apicobasal (second) and transmural (third) direction.

ldrb.calculus.bislerp(Qa: ndarray, Qb: ndarray, t: float) ndarray[source]#

Linear interpolation of two orthogonal matrices. Assume that \(Q_a\) and \(Q_b\) refers to timepoint \(0\) and \(1\) respectively. Using spherical linear interpolation (slerp) find the orthogonal matrix at timepoint \(t\).

ldrb.calculus.normalize(u: ndarray) ndarray[source]#

L2-Normalize vector with

Parameters:

u (np.ndarray) – Vector

Returns:

Normalized vector

Return type:

np.ndarray

ldrb.calculus.orient(Q: ndarray, alpha: float, beta: float) ndarray[source]#

Define the orthotropic fiber orientations.

Given a coordinate system \(Q\), in the canonical basis, rotate it in order to align with the fiber, sheet and sheet-normal axis determine by the angles \(\alpha\) (fiber) and \(\beta\) (sheets).

ldrb.calculus.quat2rot(q: ndarray) ndarray[source]#

Convert quaternion to rotation matrix

Parameters:

q (np.ndarray) – Quaternion

Returns:

Rotation matrix

Return type:

np.ndarray

ldrb.calculus.slerp(q1: ndarray, q2: ndarray, t: float) ndarray[source]#

Spherical linear interpolation from q1 to q2 at t

Parameters:
  • q1 (np.ndarray) – Source quaternion

  • q2 (np.ndarray) – Target quaternion

  • t (float) – Interpolation factor, between 0 and 1

Returns:

The spherical linear interpolation between q1 and q2 at t

Return type:

np.ndarray

class ldrb.cli.Geometry(mesh: dolfinx.mesh.Mesh, ffun: <function meshtags at 0x7f73158cf130>, markers: dict[str, list[int]])[source]#
class ldrb.io.FiberSheetSystem(f0, s0, n0)[source]#
f0: ndarray#

Alias for field number 0

n0: ndarray#

Alias for field number 2

s0: ndarray#

Alias for field number 1

class ldrb.io.LDRBOutput(f0, s0, n0, lv, rv, epi, lv_rv, apex, lv_scalar, rv_scalar, epi_scalar, lv_rv_scalar, apex_scalar, lv_gradient, rv_gradient, epi_gradient, lv_rv_gradient, apex_gradient, markers_scalar)[source]#
apex: Function | None#

Alias for field number 7

apex_gradient: Function | None#

Alias for field number 17

apex_scalar: Function | None#

Alias for field number 12

epi: Function | None#

Alias for field number 5

epi_gradient: Function | None#

Alias for field number 15

epi_scalar: Function | None#

Alias for field number 10

f0: Function#

Alias for field number 0

lv: Function | None#

Alias for field number 3

lv_gradient: Function | None#

Alias for field number 13

lv_rv: Function | None#

Alias for field number 6

lv_rv_gradient: Function | None#

Alias for field number 16

lv_rv_scalar: Function | None#

Alias for field number 11

lv_scalar: Function | None#

Alias for field number 8

markers_scalar: Function | None#

Alias for field number 18

n0: Function#

Alias for field number 2

rv: Function | None#

Alias for field number 4

rv_gradient: Function | None#

Alias for field number 14

rv_scalar: Function | None#

Alias for field number 9

s0: Function#

Alias for field number 1

ldrb.ldrb.apex_to_base(mesh: Mesh, base_marker: list[int], ffun: MeshTags) Function[source]#

Find the apex coordinate and compute the laplace equation to find the apex to base solution

Parameters:
  • mesh (dolfin.Mesh) – The mesh

  • base_marker (int) – The marker value for the basal facets

  • ffun (dolfin.MeshFunctionSizet (optional)) – A facet function containing markers for the boundaries. If not provided, the markers stored within the mesh will be used.

ldrb.ldrb.compute_fiber_sheet_system(lv_scalar: ndarray, lv_gradient: ndarray, epi_scalar: ndarray, epi_gradient: ndarray, apex_gradient: ndarray, dofs: ndarray | None = None, rv_scalar: ndarray | None = None, rv_gradient: ndarray | None = None, lv_rv_scalar: ndarray | None = None, marker_scalar: ndarray | None = None, alpha_endo_lv: float = 40, alpha_epi_lv: float = -50, alpha_endo_rv: float | None = None, alpha_epi_rv: float | None = None, alpha_endo_sept: float | None = None, alpha_epi_sept: float | None = None, beta_endo_lv: float = -65, beta_epi_lv: float = 25, beta_endo_rv: float | None = None, beta_epi_rv: float | None = None, beta_endo_sept: float | None = None, beta_epi_sept: float | None = None, epi_only: bool = False) FiberSheetSystem[source]#

Compute the fiber-sheets system on all degrees of freedom.

ldrb.ldrb.dofs_from_function_space(mesh: Mesh, fiber_space: str) ndarray[source]#

Get the dofs from a function spaces define in the fiber_space string.

ldrb.ldrb.dolfinx_ldrb(mesh: Mesh, fiber_space: str = 'P_1', ffun: MeshTags | None = None, markers: dict[str, list[int]] | dict[str, tuple[int, ...]] | dict[str, int] | None = None, alpha_endo_lv: float = 40, alpha_epi_lv: float = -50, alpha_endo_rv: float | None = None, alpha_epi_rv: float | None = None, alpha_endo_sept: float | None = None, alpha_epi_sept: float | None = None, beta_endo_lv: float = -65, beta_epi_lv: float = 25, beta_endo_rv: float | None = None, beta_epi_rv: float | None = None, beta_endo_sept: float | None = None, beta_epi_sept: float | None = None, epi_only: bool = False, **kwargs) LDRBOutput[source]#

Create fiber, cross fibers and sheet directions

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

  • fiber_space (str) – A string on the form {family}_{degree} which determines for what space the fibers should be calculated for. If not provided, then a first order Lagrange space will be used, i.e Lagrange_1.

  • ffun (dolfinx.mesh.MeshTags) – A facet function containing markers for the boundaries.

  • markers (dict[str, int | list[int]] (optional)) – A dictionary with the markers for the different boundaries defined in the facet function or within the mesh itself. The following markers must be provided: ‘base’, ‘lv’, ‘epi, ‘rv’ (optional). If the markers are not provided the following default vales will be used: base = 10, rv = 20, lv = 30, epi = 40

  • save_markers (bool) – If true save markings of the geometry. This is nice if you want to see that the LV, RV and Septum are marked correctly.

  • angles (kwargs) –

    Keyword arguments with the fiber and sheet angles. It is possible to set different angles on the LV, RV and septum, however it either the RV or septum angles are not provided, then the angles on the LV will be used. The default values are taken from the original paper, namely

    \[\begin{split}\alpha_{\text{endo}} &= 40 \\ \alpha_{\text{epi}} &= -50 \\ \beta_{\text{endo}} &= -65 \\ \beta_{\text{epi}} &= 25\end{split}\]

    The following keyword arguments are possible:

    alpha_endo_lvscalar

    Fiber angle at the LV endocardium.

    alpha_epi_lvscalar

    Fiber angle at the LV epicardium.

    beta_endo_lvscalar

    Sheet angle at the LV endocardium.

    beta_epi_lvscalar

    Sheet angle at the LV epicardium.

    alpha_endo_rvscalar

    Fiber angle at the RV endocardium.

    alpha_epi_rvscalar

    Fiber angle at the RV epicardium.

    beta_endo_rvscalar

    Sheet angle at the RV endocardium.

    beta_epi_rvscalar

    Sheet angle at the RV epicardium.

    alpha_endo_septscalar

    Fiber angle at the septum endocardium.

    alpha_epi_septscalar

    Fiber angle at the septum epicardium.

    beta_endo_septscalar

    Sheet angle at the septum endocardium.

    beta_epi_septscalar

    Sheet angle at the septum epicardium.

  • epi_only (bool) – If true, then only the angles on the epicardial surface will be used. This is a hack to make it possible to compute purely longitudinal fibers

ldrb.ldrb.project_gradients(mesh: Mesh, scalar_solutions: dict[str, Function], fiber_space: str = 'P_1') tuple[dict[str, ndarray], dict[str, Function]][source]#

Calculate the gradients using projections

Parameters:
  • mesh (dolfin.Mesh) – The mesh

  • fiber_space (str) – A string on the form {family}_{degree} which determines for what space the fibers should be calculated for.

  • scalar_solutions (dict) – A dictionary with the scalar solutions that you want to compute the gradients of.

ldrb.ldrb.scalar_laplacians(mesh: Mesh, markers: dict[str, list[int]], ffun: MeshTags) dict[str, Function][source]#

Calculate the laplacians

Parameters:
  • mesh (dolfin.Mesh) – A dolfin mesh

  • markers (dict (optional)) – A dictionary with the markers for the different boundaries defined in the facet function or within the mesh itself. The following markers must be provided: ‘base’, ‘lv’, ‘epi, ‘rv’ (optional). If the markers are not provided the following default vales will be used: base = 10, rv = 20, lv = 30, epi = 40.

ldrb.ldrb.standard_dofs(n: int) ndarray[source]#

Get the standard list of dofs for a given length

ldrb.ldrb.transform_markers(markers: dict[str, list[int]]) dict[str, list[int]][source]#

Convert markers generated with gmsh to format for ldrb

ldrb.utils.default_markers() dict[str, list[int]][source]#

Default markers for the mesh boundaries

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

Parse a string representation of a basix element family

ldrb.utils.space_from_string(space_string: str, mesh: Mesh, dim: int) 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.