ldrb package

Submodules

ldrb.calculus module

ldrb.calculus.axis(u: numpy.ndarray, v: numpy.ndarray) numpy.ndarray

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: numpy.ndarray, Qb: numpy.ndarray, t: float) numpy.ndarray

Linear interpolation of two orthogonal matrices. Assiume 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: numpy.ndarray) numpy.ndarray

Normalize vector

ldrb.calculus.orient(Q: numpy.ndarray, alpha: float, beta: float) numpy.ndarray

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.system_at_dof(lv: float, rv: float, epi: float, grad_lv: numpy.ndarray, grad_rv: numpy.ndarray, grad_epi: numpy.ndarray, grad_ab: numpy.ndarray, alpha_endo: float, alpha_epi: float, beta_endo: float, beta_epi: float, tol: float = 1e-07) numpy.ndarray

Compte the fiber, sheet and sheet normal at a single degre of freedom

Parameters
  • lv (float) – Value of the Laplace solution for the LV at the dof.

  • rv (float) – Value of the Laplace solution for the RV at the dof.

  • epi (float) – Value of the Laplace solution for the EPI at the dof.

  • grad_lv (np.ndarray) – Gradient of the Laplace solution for the LV at the dof.

  • grad_rv (np.ndarray) – Gradient of the Laplace solution for the RV at the dof.

  • grad_epi (np.ndarray) – Gradient of the Laplace solution for the EPI at the dof.

  • grad_epi – Gradient of the Laplace solution for the apex to base. at the dof

  • alpha_endo (scalar) – Fiber angle at the endocardium.

  • alpha_epi (scalar) – Fiber angle at the epicardium.

  • beta_endo (scalar) – Sheet angle at the endocardium.

  • beta_epi (scalar) – Sheet angle at the epicardium.

  • tol (scalar) – Tolerance for whether to consider the scalar values. Default: 1e-7

ldrb.ldrb module

class ldrb.ldrb.FiberSheetSystem(fiber, sheet, sheet_normal)

Bases: tuple

fiber

Alias for field number 0

sheet

Alias for field number 1

sheet_normal

Alias for field number 2

ldrb.ldrb.apex_to_base(mesh: dolfin.cpp.mesh.Mesh, base_marker: int, ffun: dolfin.mesh.meshfunction.MeshFunction, use_krylov_solver: bool = False, krylov_solver_atol: Optional[float] = None, krylov_solver_rtol: Optional[float] = None, krylov_solver_max_its: Optional[int] = None, verbose: bool = False)

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.

  • use_krylov_solver (bool) – If True use Krylov solver, by default False

  • krylov_solver_atol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the absolute residual. Default: 1e-15.

  • krylov_solver_rtol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the relative residual. Default: 1e-10.

  • krylov_solver_max_its (int (optional)) – If a Krylov solver is used, this option specifies the maximum number of iterations to perform. Default: 10000.

  • verbose (bool) – If true, print more info, by default False

ldrb.ldrb.bayer(cases, mesh, markers, ffun, verbose: bool, use_krylov_solver: bool, strict: bool, krylov_solver_atol: Optional[float] = None, krylov_solver_rtol: Optional[float] = None, krylov_solver_max_its: Optional[int] = None) Dict[str, dolfin.function.function.Function]
ldrb.ldrb.check_boundaries_are_marked(mesh: dolfin.cpp.mesh.Mesh, ffun: dolfin.mesh.meshfunction.MeshFunction, markers: Dict[str, int], boundaries: List[str]) None
ldrb.ldrb.compute_fiber_sheet_system(lv_scalar: numpy.ndarray, lv_gradient: numpy.ndarray, epi_scalar: numpy.ndarray, epi_gradient: numpy.ndarray, apex_gradient: numpy.ndarray, dofs: Optional[numpy.ndarray] = None, rv_scalar: Optional[numpy.ndarray] = None, rv_gradient: Optional[numpy.ndarray] = None, lv_rv_scalar: Optional[numpy.ndarray] = None, marker_scalar: Optional[numpy.ndarray] = None, alpha_endo_lv: float = 40, alpha_epi_lv: float = - 50, alpha_endo_rv: Optional[float] = None, alpha_epi_rv: Optional[float] = None, alpha_endo_sept: Optional[float] = None, alpha_epi_sept: Optional[float] = None, beta_endo_lv: float = - 65, beta_epi_lv: float = 25, beta_endo_rv: Optional[float] = None, beta_epi_rv: Optional[float] = None, beta_endo_sept: Optional[float] = None, beta_epi_sept: Optional[float] = None) ldrb.ldrb.FiberSheetSystem

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

ldrb.ldrb.dofs_from_function_space(mesh: dolfin.cpp.mesh.Mesh, fiber_space: str) numpy.ndarray

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

ldrb.ldrb.dolfin_ldrb(mesh: dolfin.cpp.mesh.Mesh, fiber_space: str = 'CG_1', ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, markers: Optional[Dict[str, int]] = None, log_level: int = 20, use_krylov_solver: bool = False, krylov_solver_atol: Optional[float] = None, krylov_solver_rtol: Optional[float] = None, krylov_solver_max_its: Optional[int] = None, strict: bool = False, save_markers: bool = False, alpha_endo_lv: float = 40, alpha_epi_lv: float = - 50, alpha_endo_rv: Optional[float] = None, alpha_epi_rv: Optional[float] = None, alpha_endo_sept: Optional[float] = None, alpha_epi_sept: Optional[float] = None, beta_endo_lv: float = - 65, beta_epi_lv: float = 25, beta_endo_rv: Optional[float] = None, beta_epi_rv: Optional[float] = None, beta_endo_sept: Optional[float] = None, beta_epi_sept: Optional[float] = None)

Create fiber, cross fibers and sheet directions

Parameters
  • mesh (dolfin.Mesh) – The mesh

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

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

  • markers (dict (optional)) – A dictionary with the markers for the different bondaries defined in the facet function or within the mesh itself. The follwing 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

  • log_level (int) – How much to print. DEBUG=10, INFO=20, WARNING=30. Default: INFO

  • use_krylov_solver (bool) – If True use Krylov solver, by default False

  • krylov_solver_atol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the absolute residual. Default: 1e-15.

  • krylov_solver_rtol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the relative residual. Default: 1e-10.

  • krylov_solver_max_its (int (optional)) – If a Krylov solver is used, this option specifies the maximum number of iterations to perform. Default: 10000.

  • strict (bool) – If true raise RuntimeError if solutions does not sum to 1.0

  • 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.

ldrb.ldrb.fiber_system_to_dolfin(system: ldrb.ldrb.FiberSheetSystem, mesh: dolfin.cpp.mesh.Mesh, fiber_space: str) ldrb.ldrb.FiberSheetSystem

Convert fiber-sheet system of numpy arrays to dolfin functions.

ldrb.ldrb.find_cases_and_boundaries(ffun: dolfin.mesh.meshfunction.MeshFunction, markers: Optional[Dict[str, int]]) Tuple[List[str], List[str], Dict[str, int]]
ldrb.ldrb.laplace(mesh: dolfin.cpp.mesh.Mesh, markers: Optional[Dict[str, int]], fiber_space: str = 'CG_1', ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, use_krylov_solver: bool = False, krylov_solver_atol: Optional[float] = None, krylov_solver_rtol: Optional[float] = None, krylov_solver_max_its: Optional[int] = None, verbose: bool = False, strict: bool = False) Dict[str, numpy.ndarray]

Solve the laplace equation and project the gradients of the solutions.

ldrb.ldrb.project_gradients(mesh: dolfin.cpp.mesh.Mesh, scalar_solutions: Dict[str, dolfin.function.function.Function], fiber_space: str = 'CG_1') Dict[str, numpy.ndarray]

Calculate the gradients using projections

Parameters
  • mesh (dolfin.Mesh) – The mesh

  • fiber_space (str) – A string on the form {familiy}_{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: dolfin.cpp.mesh.Mesh, markers: Optional[Dict[str, int]] = None, ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, use_krylov_solver: bool = False, krylov_solver_atol: Optional[float] = None, krylov_solver_rtol: Optional[float] = None, krylov_solver_max_its: Optional[int] = None, verbose: bool = False, strict: bool = False) Dict[str, dolfin.function.function.Function]

Calculate the laplacians

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

  • markers (dict (optional)) – A dictionary with the markers for the different bondaries defined in the facet function or within the mesh itself. The follwing 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.

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

  • use_krylov_solver (bool) – If True use Krylov solver, by default False

  • krylov_solver_atol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the absolute residual. Default: 1e-15.

  • krylov_solver_rtol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the relative residual. Default: 1e-10.

  • krylov_solver_max_its (int (optional)) – If a Krylov solver is used, this option specifies the maximum number of iterations to perform. Default: 10000.

  • verbose (bool) – If true, print more info, by default False

  • strict (bool) – If true raise RuntimeError if solutions does not sum to 1.0

ldrb.ldrb.solve_krylov(a, L, bcs, u: dolfin.function.function.Function, verbose: bool = False, ksp_type='cg', ksp_norm_type='unpreconditioned', ksp_atol=1e-15, ksp_rtol=1e-10, ksp_max_it=10000, ksp_error_if_not_converged=False, pc_type='hypre') dolfin.cpp.la.PETScKrylovSolver
ldrb.ldrb.solve_regular(a, L, bcs, u, solver_parameters)
ldrb.ldrb.solve_system(a, L, bcs, u: dolfin.function.function.Function, solver_parameters: Optional[Dict[str, str]] = None, use_krylov_solver: bool = False, krylov_solver_atol: Optional[float] = None, krylov_solver_rtol: Optional[float] = None, krylov_solver_max_its: Optional[int] = None, verbose: bool = False) Optional[dolfin.cpp.la.PETScKrylovSolver]
ldrb.ldrb.standard_dofs(n: int) numpy.ndarray

Get the standard list of dofs for a given length

ldrb.save module

ldrb.save.dolfin_to_hd5(obj, h5name, time='', comm=<mpi4py.MPI.Intracomm object>, name=None)

Save object to and HDF file.

Parameters
  • obj (dolfin.Mesh or dolfin.Function) – The object you want to save

  • name (str) – Name of the object

  • h5group (str) – The folder you want to save the object withing the HDF file. Default: ‘’

ldrb.save.fiber_to_xdmf(fun, fname, comm=<mpi4py.MPI.Intracomm object>)
ldrb.save.fun_to_xdmf(fun, fname, name='function')
ldrb.save.load_dict_from_h5(fname, h5group='', comm=<mpi4py.MPI.Intracomm object>)

Load the given h5file into a dictionary

ldrb.save.save_scalar_function(comm, obj, h5name, h5group='', file_mode='w')
ldrb.save.save_vector_function(comm, obj, h5name, h5group='', file_mode='w')

ldrb.utils module

class ldrb.utils.Geometry(mesh, ffun, markers)

Bases: tuple

ffun

Alias for field number 1

markers

Alias for field number 2

mesh

Alias for field number 0

ldrb.utils.assign_to_vector(v, a)

Assign the value of the array a to the dolfin vector v

ldrb.utils.broadcast(array, from_process)

Broadcast array to all processes

ldrb.utils.convert_msh_to_xdmf(msh_file, triangle_mesh_name, tetra_mesh_name)
ldrb.utils.create_biv_mesh(N: int = 13, a_endo_lv: float = 1.5, b_endo_lv: float = 0.5, c_endo_lv: float = 0.5, a_epi_lv: float = 2.0, b_epi_lv: float = 1.0, c_epi_lv: float = 1.0, center_lv: Tuple[float, float, float] = (0.0, 0.0, 0.0), a_endo_rv: float = 1.45, b_endo_rv: float = 1.25, c_endo_rv: float = 0.75, a_epi_rv: float = 1.75, b_epi_rv: float = 1.5, c_epi_rv: float = 1.0, center_rv: Tuple[float, float, float] = (0.0, 0.5, 0.0), base_x: float = 0.0, markers: Optional[Dict[str, int]] = None) ldrb.utils.Geometry

Create an biv-ellipsoidal mesh.

An ellipsoid is given by the equation

\[\frac{x^2}{a} + \frac{y^2}{b} + \frac{z^2}{c} = 1\]

We create three ellipsoids, one for the LV and RV endocardium and one for the epicardium and subtract them and then cut the base. For simplicity we assume that the longitudinal axis is in in \(x\)-direction and as default the base is located at the \(x=0\) plane.

ldrb.utils.create_lv_mesh(N: int = 13, a_endo: float = 1.5, b_endo: float = 0.5, c_endo: float = 0.5, a_epi: float = 2.0, b_epi: float = 1.0, c_epi: float = 1.0, center: Tuple[float, float, float] = (0.0, 0.0, 0.0), base_x: float = 0.0, markers: Optional[Dict[str, int]] = None) ldrb.utils.Geometry

Create an lv-ellipsoidal mesh.

An ellipsoid is given by the equation

\[\frac{x^2}{a} + \frac{y^2}{b} + \frac{z^2}{c} = 1\]

We create two ellipsoids, one for the endocardium and one for the epicardium and subtract them and then cut the base. For simplicity we assume that the longitudinal axis is in in \(x\)-direction and as default the base is located at the \(x=0\) plane.

ldrb.utils.create_mesh(mesh, cell_type)
ldrb.utils.default_markers() Dict[str, int]

Default markers for the mesh boundaries

ldrb.utils.distribution(number)

Get distribution of number on all processes

ldrb.utils.gather(array, on_process=0, flatten=False)

Gather array from all processes on a single process

ldrb.utils.gather_broadcast(arr)
ldrb.utils.gmsh2dolfin(msh_file, unlink: bool = True)
ldrb.utils.has_meshio() bool
ldrb.utils.has_mshr()
ldrb.utils.mark_biv_mesh(mesh: dolfin.cpp.mesh.Mesh, ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, markers: Optional[Dict[str, int]] = None, tol: float = 0.01, values: Dict[str, int] = {'lv': 0, 'rv': 2, 'septum': 1}) dolfin.mesh.meshfunction.MeshFunction
ldrb.utils.mark_facets(mesh: dolfin.cpp.mesh.Mesh, ffun: dolfin.mesh.meshfunction.MeshFunction)

Mark mesh according to facet function

ldrb.utils.mpi_comm_world()
ldrb.utils.read_meshfunction(fname, obj)
ldrb.utils.space_from_string(space_string: str, mesh: dolfin.cpp.mesh.Mesh, dim: int) dolfin.function.functionspace.FunctionSpace

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

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

  • mesh (dolfin.Mesh) – The mesh

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

ldrb.utils.value_size(obj: ufl.coefficient.Coefficient) Union[List[int], int]

Module contents

ldrb.create_biv_mesh(N: int = 13, a_endo_lv: float = 1.5, b_endo_lv: float = 0.5, c_endo_lv: float = 0.5, a_epi_lv: float = 2.0, b_epi_lv: float = 1.0, c_epi_lv: float = 1.0, center_lv: Tuple[float, float, float] = (0.0, 0.0, 0.0), a_endo_rv: float = 1.45, b_endo_rv: float = 1.25, c_endo_rv: float = 0.75, a_epi_rv: float = 1.75, b_epi_rv: float = 1.5, c_epi_rv: float = 1.0, center_rv: Tuple[float, float, float] = (0.0, 0.5, 0.0), base_x: float = 0.0, markers: Optional[Dict[str, int]] = None) ldrb.utils.Geometry

Create an biv-ellipsoidal mesh.

An ellipsoid is given by the equation

\[\frac{x^2}{a} + \frac{y^2}{b} + \frac{z^2}{c} = 1\]

We create three ellipsoids, one for the LV and RV endocardium and one for the epicardium and subtract them and then cut the base. For simplicity we assume that the longitudinal axis is in in \(x\)-direction and as default the base is located at the \(x=0\) plane.

ldrb.create_lv_mesh(N: int = 13, a_endo: float = 1.5, b_endo: float = 0.5, c_endo: float = 0.5, a_epi: float = 2.0, b_epi: float = 1.0, c_epi: float = 1.0, center: Tuple[float, float, float] = (0.0, 0.0, 0.0), base_x: float = 0.0, markers: Optional[Dict[str, int]] = None) ldrb.utils.Geometry

Create an lv-ellipsoidal mesh.

An ellipsoid is given by the equation

\[\frac{x^2}{a} + \frac{y^2}{b} + \frac{z^2}{c} = 1\]

We create two ellipsoids, one for the endocardium and one for the epicardium and subtract them and then cut the base. For simplicity we assume that the longitudinal axis is in in \(x\)-direction and as default the base is located at the \(x=0\) plane.

ldrb.dolfin_ldrb(mesh: dolfin.cpp.mesh.Mesh, fiber_space: str = 'CG_1', ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, markers: Optional[Dict[str, int]] = None, log_level: int = 20, use_krylov_solver: bool = False, krylov_solver_atol: Optional[float] = None, krylov_solver_rtol: Optional[float] = None, krylov_solver_max_its: Optional[int] = None, strict: bool = False, save_markers: bool = False, alpha_endo_lv: float = 40, alpha_epi_lv: float = - 50, alpha_endo_rv: Optional[float] = None, alpha_epi_rv: Optional[float] = None, alpha_endo_sept: Optional[float] = None, alpha_epi_sept: Optional[float] = None, beta_endo_lv: float = - 65, beta_epi_lv: float = 25, beta_endo_rv: Optional[float] = None, beta_epi_rv: Optional[float] = None, beta_endo_sept: Optional[float] = None, beta_epi_sept: Optional[float] = None)

Create fiber, cross fibers and sheet directions

Parameters
  • mesh (dolfin.Mesh) – The mesh

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

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

  • markers (dict (optional)) – A dictionary with the markers for the different bondaries defined in the facet function or within the mesh itself. The follwing 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

  • log_level (int) – How much to print. DEBUG=10, INFO=20, WARNING=30. Default: INFO

  • use_krylov_solver (bool) – If True use Krylov solver, by default False

  • krylov_solver_atol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the absolute residual. Default: 1e-15.

  • krylov_solver_rtol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the relative residual. Default: 1e-10.

  • krylov_solver_max_its (int (optional)) – If a Krylov solver is used, this option specifies the maximum number of iterations to perform. Default: 10000.

  • strict (bool) – If true raise RuntimeError if solutions does not sum to 1.0

  • 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.

ldrb.fiber_to_xdmf(fun, fname, comm=<mpi4py.MPI.Intracomm object>)
ldrb.fun_to_xdmf(fun, fname, name='function')
ldrb.gmsh2dolfin(msh_file, unlink: bool = True)
ldrb.project_gradients(mesh: dolfin.cpp.mesh.Mesh, scalar_solutions: Dict[str, dolfin.function.function.Function], fiber_space: str = 'CG_1') Dict[str, numpy.ndarray]

Calculate the gradients using projections

Parameters
  • mesh (dolfin.Mesh) – The mesh

  • fiber_space (str) – A string on the form {familiy}_{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: dolfin.cpp.mesh.Mesh, markers: Optional[Dict[str, int]] = None, ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, use_krylov_solver: bool = False, krylov_solver_atol: Optional[float] = None, krylov_solver_rtol: Optional[float] = None, krylov_solver_max_its: Optional[int] = None, verbose: bool = False, strict: bool = False) Dict[str, dolfin.function.function.Function]

Calculate the laplacians

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

  • markers (dict (optional)) – A dictionary with the markers for the different bondaries defined in the facet function or within the mesh itself. The follwing 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.

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

  • use_krylov_solver (bool) – If True use Krylov solver, by default False

  • krylov_solver_atol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the absolute residual. Default: 1e-15.

  • krylov_solver_rtol (float (optional)) – If a Krylov solver is used, this option specifies a convergence criterion in terms of the relative residual. Default: 1e-10.

  • krylov_solver_max_its (int (optional)) – If a Krylov solver is used, this option specifies the maximum number of iterations to perform. Default: 10000.

  • verbose (bool) – If true, print more info, by default False

  • strict (bool) – If true raise RuntimeError if solutions does not sum to 1.0

ldrb.space_from_string(space_string: str, mesh: dolfin.cpp.mesh.Mesh, dim: int) dolfin.function.functionspace.FunctionSpace

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

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

  • mesh (dolfin.Mesh) – The mesh

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