ldrb package#
Submodules#
ldrb.calculus module#
- ldrb.calculus.axis(u: ndarray, v: ndarray) 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: ndarray, Qb: ndarray, t: float) ndarray #
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.compute_fiber_sheet_system(f0, s0, n0, xdofs, ydofs, zdofs, sdofs, lv_scalar, rv_scalar, epi_scalar, lv_rv_scalar, lv_gradient, rv_gradient, epi_gradient, apex_gradient, marker_scalar, alpha_endo_lv, alpha_epi_lv, alpha_endo_rv, alpha_epi_rv, alpha_endo_sept, alpha_epi_sept, beta_endo_lv, beta_epi_lv, beta_endo_rv, beta_epi_rv, beta_endo_sept, beta_epi_sept)#
- ldrb.calculus.normalize(u: ndarray) ndarray #
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 #
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 #
Convert quaternion to rotation matrix
- Parameters:
q (np.ndarray) – Quaternion
- Returns:
Rotation matrix
- Return type:
np.ndarray
- ldrb.calculus.rot2quat(Q: ndarray) ndarray #
- ldrb.calculus.slerp(q1: ndarray, q2: ndarray, t: float) ndarray #
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
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: Mesh, base_marker: List[int], ffun: MeshFunction, solver: PETScKrylovSolver) Function #
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.
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: List[str], mesh: Mesh, markers: Dict[str, List[int]], ffun: MeshFunction, verbose: bool, strict: bool, krylov_solver_atol: float | None = None, krylov_solver_rtol: float | None = None, krylov_solver_max_its: int | None = None, ksp_monitor: bool = False, ksp_view: bool = False, pc_view: bool = False) Dict[str, Function] #
- ldrb.ldrb.check_boundaries_are_marked(mesh: Mesh, ffun: MeshFunction, markers: Dict[str, List[int]], boundaries: List[str]) None #
- 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) FiberSheetSystem #
Compute the fiber-sheets system on all degrees of freedom.
- ldrb.ldrb.dofs_from_function_space(mesh: Mesh, fiber_space: str) ndarray #
Get the dofs from a function spaces define in the fiber_space string.
- ldrb.ldrb.dolfin_ldrb(mesh: Mesh, fiber_space: str = 'CG_1', ffun: MeshFunction | None = None, markers: Dict[str, int | List[int]] | None = None, log_level: int = 20, krylov_solver_atol: float | None = None, krylov_solver_rtol: float | None = None, krylov_solver_max_its: int | None = None, strict: bool = False, save_markers: bool = False, 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) FiberSheetSystem #
Create fiber, cross fibers and sheet directions
- 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. If not provided, 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[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
log_level (int) – How much to print. DEBUG=10, INFO=20, WARNING=30. Default: INFO
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: FiberSheetSystem, mesh: Mesh, fiber_space: str) FiberSheetSystem #
Convert fiber-sheet system of numpy arrays to dolfin functions.
- ldrb.ldrb.find_cases_and_boundaries(markers_to_process: Dict[str, int | List[int]] | None) Tuple[List[str], List[str], Dict[str, List[int]]] #
- ldrb.ldrb.laplace(mesh: Mesh, markers_to_process: Dict[str, int | List[int]] | None, fiber_space: str = 'CG_1', ffun: MeshFunction | None = None, krylov_solver_atol: float | None = None, krylov_solver_rtol: float | None = None, krylov_solver_max_its: int | None = None, verbose: bool = False, strict: bool = False) Dict[str, ndarray] #
Solve the laplace equation and project the gradients of the solutions.
- ldrb.ldrb.process_markers(markers: Dict[str, int | List[int]] | None) Dict[str, List[int]] #
- ldrb.ldrb.project_gradients(mesh: Mesh, scalar_solutions: Dict[str, Function], fiber_space: str = 'CG_1') Dict[str, ndarray] #
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_to_process: dict[str, int | list[int]] | None = None, ffun: MeshFunction | None = None, krylov_solver_atol: float | None = None, krylov_solver_rtol: float | None = None, krylov_solver_max_its: int | None = None, verbose: bool = False, strict: bool = False) Dict[str, Function] #
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.
fiber_space (str) – A string on the form {family}_{degree} which determines for what space the fibers should be calculated for.
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.standard_dofs(n: int) 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.Projector(V: FunctionSpace, solver_type: str = 'lu', preconditioner_type: str = 'default')#
Bases:
object
- project(u: Function, f: Expr) None #
Project f into u. :param u: The function to project into :type u: df.Function :param f: The ufl expression to project :type f: ufl.core.expr.Expr
- 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.default_markers() Dict[str, List[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.mpi_comm_world()#
- ldrb.utils.space_from_string(space_string: str, mesh: Mesh, dim: int) FunctionSpace #
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.utils.value_size(obj: Coefficient) List[int] | int #
Module contents#
- ldrb.dolfin_ldrb(mesh: Mesh, fiber_space: str = 'CG_1', ffun: MeshFunction | None = None, markers: Dict[str, int | List[int]] | None = None, log_level: int = 20, krylov_solver_atol: float | None = None, krylov_solver_rtol: float | None = None, krylov_solver_max_its: int | None = None, strict: bool = False, save_markers: bool = False, 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) FiberSheetSystem #
Create fiber, cross fibers and sheet directions
- 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. If not provided, 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[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
log_level (int) – How much to print. DEBUG=10, INFO=20, WARNING=30. Default: INFO
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.project_gradients(mesh: Mesh, scalar_solutions: Dict[str, Function], fiber_space: str = 'CG_1') Dict[str, ndarray] #
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_to_process: dict[str, int | list[int]] | None = None, ffun: MeshFunction | None = None, krylov_solver_atol: float | None = None, krylov_solver_rtol: float | None = None, krylov_solver_max_its: int | None = None, verbose: bool = False, strict: bool = False) Dict[str, Function] #
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.
fiber_space (str) – A string on the form {family}_{degree} which determines for what space the fibers should be calculated for.
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: Mesh, dim: int) FunctionSpace #
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.