API Reference#

class cardiac_geometries.Geometry(mesh: dolfinx.mesh.Mesh, markers: dict[str, tuple[int, int]] = <factory>, cfun: dolfinx.mesh.MeshTags | None = None, ffun: dolfinx.mesh.MeshTags | None = None, efun: dolfinx.mesh.MeshTags | None = None, vfun: dolfinx.mesh.MeshTags | None = None, f0: dolfinx.fem.function.Function | None = None, s0: dolfinx.fem.function.Function | None = None, n0: dolfinx.fem.function.Function | None = None, info: dict[str, typing.Any] = <factory>)[source]#
property ds#

Surface measure for the mesh using the facet function ffun if it exists as subdomain data.

property dx#

Volume measure for the mesh using the cell function cfun if it exists as subdomain data.

property facet_normal: FacetNormal#

Facet normal vector for the mesh.

classmethod from_file(comm: Intracomm, path: str | Path, function_space_data: dict[str, ndarray] | None = None) Geometry[source]#

Read geometry from a file using adios4dolfinx.

Parameters:
  • comm (MPI.Intracomm) – The MPI communicator to use for reading the mesh.

  • path (str | Path) – The path to the file containing the geometry data.

  • function_space_data (dict[str, np.ndarray] | None, optional) – A dictionary containing function space data for the functions to be read. If None, it will be read from the file.

Returns:

An instance of the Geometry class containing the mesh, markers, and functions.

Return type:

Geometry

classmethod from_folder(comm: Intracomm, folder: str | Path) Geometry[source]#

Read geometry from a folder containing mesh and markers files.

Parameters:
  • comm (MPI.Intracomm) – The MPI communicator to use for reading the mesh and markers.

  • folder (str | Path) – The path to the folder containing the geometry data. The folder should contain the following files: - mesh.xdmf: The mesh file in XDMF format. - markers.json: A JSON file containing markers. - microstructure.json: A JSON file containing microstructure data (optional). - microstructure.bp: A BP file containing microstructure functions (optional). - info.json: A JSON file containing additional information (optional).

Returns:

An instance of the Geometry class containing the mesh, markers, and functions.

Return type:

Geometry

Raises:

ValueError – If the required mesh file is not found in the specified folder.

refine(n=1, outdir: Path | None = None) Geometry[source]#

Refine the mesh and transfer the meshtags to new geometry. Also regenerate fibers if self.info is found. If self.info is not found, it currently raises a NotImplementedError, however fiber could be interpolated from the old mesh to the new mesh but this will result in a loss of information about the fiber orientation.

Parameters:
  • n (int, optional) – Number of times to refine the mesh, by default 1

  • outdir (Path | None, optional) – Output directory to save the refined mesh and meshtags, by default None in which case the mesh is not saved.

Returns:

A new Geometry object with the refined mesh and updated meshtags.

Return type:

Geometry

Raises:

NotImplementedError – If self.info is not found, indicating that fiber interpolation after refinement is not implemented yet.

save(path: str | Path) None[source]#

Save the geometry to a file using adios4dolfinx.

Parameters:

path (str | Path) – The path to the file where the geometry will be saved. The file will be created if it does not exist, or overwritten if it does.

class cardiac_geometries.geometry.Geometry(mesh: dolfinx.mesh.Mesh, markers: dict[str, tuple[int, int]] = <factory>, cfun: dolfinx.mesh.MeshTags | None = None, ffun: dolfinx.mesh.MeshTags | None = None, efun: dolfinx.mesh.MeshTags | None = None, vfun: dolfinx.mesh.MeshTags | None = None, f0: dolfinx.fem.function.Function | None = None, s0: dolfinx.fem.function.Function | None = None, n0: dolfinx.fem.function.Function | None = None, info: dict[str, typing.Any] = <factory>)[source]#
property ds#

Surface measure for the mesh using the facet function ffun if it exists as subdomain data.

property dx#

Volume measure for the mesh using the cell function cfun if it exists as subdomain data.

property facet_normal: FacetNormal#

Facet normal vector for the mesh.

classmethod from_file(comm: Intracomm, path: str | Path, function_space_data: dict[str, ndarray] | None = None) Geometry[source]#

Read geometry from a file using adios4dolfinx.

Parameters:
  • comm (MPI.Intracomm) – The MPI communicator to use for reading the mesh.

  • path (str | Path) – The path to the file containing the geometry data.

  • function_space_data (dict[str, np.ndarray] | None, optional) – A dictionary containing function space data for the functions to be read. If None, it will be read from the file.

Returns:

An instance of the Geometry class containing the mesh, markers, and functions.

Return type:

Geometry

classmethod from_folder(comm: Intracomm, folder: str | Path) Geometry[source]#

Read geometry from a folder containing mesh and markers files.

Parameters:
  • comm (MPI.Intracomm) – The MPI communicator to use for reading the mesh and markers.

  • folder (str | Path) – The path to the folder containing the geometry data. The folder should contain the following files: - mesh.xdmf: The mesh file in XDMF format. - markers.json: A JSON file containing markers. - microstructure.json: A JSON file containing microstructure data (optional). - microstructure.bp: A BP file containing microstructure functions (optional). - info.json: A JSON file containing additional information (optional).

Returns:

An instance of the Geometry class containing the mesh, markers, and functions.

Return type:

Geometry

Raises:

ValueError – If the required mesh file is not found in the specified folder.

refine(n=1, outdir: Path | None = None) Geometry[source]#

Refine the mesh and transfer the meshtags to new geometry. Also regenerate fibers if self.info is found. If self.info is not found, it currently raises a NotImplementedError, however fiber could be interpolated from the old mesh to the new mesh but this will result in a loss of information about the fiber orientation.

Parameters:
  • n (int, optional) – Number of times to refine the mesh, by default 1

  • outdir (Path | None, optional) – Output directory to save the refined mesh and meshtags, by default None in which case the mesh is not saved.

Returns:

A new Geometry object with the refined mesh and updated meshtags.

Return type:

Geometry

Raises:

NotImplementedError – If self.info is not found, indicating that fiber interpolation after refinement is not implemented yet.

save(path: str | Path) None[source]#

Save the geometry to a file using adios4dolfinx.

Parameters:

path (str | Path) – The path to the file where the geometry will be saved. The file will be created if it does not exist, or overwritten if it does.

cardiac_geometries.mesh.biv_ellipsoid(outdir: str | ~pathlib.Path, char_length: float = 0.5, center_lv_x: float = 0.0, center_lv_y: float = 0.0, center_lv_z: float = 0.0, a_endo_lv: float = 2.5, b_endo_lv: float = 1.0, c_endo_lv: float = 1.0, a_epi_lv: float = 3.0, b_epi_lv: float = 1.5, c_epi_lv: float = 1.5, center_rv_x: float = 0.0, center_rv_y: float = 0.5, center_rv_z: float = 0.0, a_endo_rv: float = 3.0, b_endo_rv: float = 1.5, c_endo_rv: float = 1.5, a_epi_rv: float = 4.0, b_epi_rv: float = 2.5, c_epi_rv: float = 2.0, create_fibers: bool = False, fiber_angle_endo: float = 60, fiber_angle_epi: float = -60, fiber_space: str = 'P_1', verbose: bool = False, comm: ~mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>) Geometry[source]#

Create BiV ellipsoidal geometry

Parameters:
  • outdir (str | Path) – Directory where to save the results.

  • char_length (float, optional) – Characteristic length of mesh, by default 0.5

  • center_lv_y (float, optional) – X-coordinate for the center of the lv, by default 0.0

  • center_lv_y – Y-coordinate for the center of the lv, by default 0.0

  • center_lv_z (float, optional) – Z-coordinate for the center of the lv, by default 0.0

  • a_endo_lv (float, optional) – Dilation of lv endo ellipsoid in the x-direction, by default 2.5

  • b_endo_lv (float, optional) – Dilation of lv endo ellipsoid in the y-direction, by default 1.0

  • c_endo_lv (float, optional) – Dilation of lv endo ellipsoid in the z-direction, by default 1.0

  • a_epi_lv (float, optional) – Dilation of lv epi ellipsoid in the x-direction, by default 3.0

  • b_epi_lv (float, optional) – Dilation of lv epi ellipsoid in the y-direction, by default 1.5

  • c_epi_lv (float, optional) – Dilation of lv epi ellipsoid in the z-direction, by default 1.5

  • center_rv_x (float, optional) – X-coordinate for the center of the rv, by default 0.0

  • center_rv_y (float, optional) – Y-coordinate for the center of the rv, by default 0.5

  • center_rv_z (float, optional) – Z-coordinate for the center of the rv, by default 0.0

  • a_endo_rv (float, optional) – Dilation of rv endo ellipsoid in the x-direction, by default 3.0

  • b_endo_rv (float, optional) – Dilation of rv endo ellipsoid in the y-direction, by default 1.5

  • c_endo_rv (float, optional) – Dilation of rv endo ellipsoid in the z-direction, by default 1.5

  • a_epi_rv (float, optional) – Dilation of rv epi ellipsoid in the x-direction, by default 4.0

  • b_epi_rv (float, optional) – Dilation of rv epi ellipsoid in the y-direction, by default 2.5

  • c_epi_rv (float, optional) – Dilation of rv epi ellipsoid in the z-direction, by default 2.0

  • create_fibers (bool, optional) – If True create analytic fibers, by default False

  • fiber_angle_endo (float, optional) – Angle for the endocardium, by default 60

  • fiber_angle_epi (float, optional) – Angle for the epicardium, by default -60

  • fiber_space (str, optional) – Function space for fibers of the form family_degree, by default “P_1”

  • verbose (bool, optional) – If True print information from gmsh, by default False

  • comm (MPI.Comm, optional) – MPI communicator, by default MPI.COMM_WORLD

Returns:

A Geometry with the mesh, markers, markers functions and fibers.

Return type:

cardiac_geometries.geometry.Geometry

cardiac_geometries.mesh.biv_ellipsoid_torso(outdir: str | ~pathlib.Path, char_length: float = 0.5, heart_as_surface: bool = False, torso_length: float = 20.0, torso_width: float = 20.0, torso_height: float = 20.0, rotation_angle: float = 0.5235987755982988, center_lv_x: float = 0.0, center_lv_y: float = 0.0, center_lv_z: float = 0.0, a_endo_lv: float = 2.5, b_endo_lv: float = 1.0, c_endo_lv: float = 1.0, a_epi_lv: float = 3.0, b_epi_lv: float = 1.5, c_epi_lv: float = 1.5, center_rv_x: float = 0.0, center_rv_y: float = 0.5, center_rv_z: float = 0.0, a_endo_rv: float = 3.0, b_endo_rv: float = 1.5, c_endo_rv: float = 1.5, a_epi_rv: float = 4.0, b_epi_rv: float = 2.5, c_epi_rv: float = 2.0, create_fibers: bool = False, fiber_angle_endo: float = 60, fiber_angle_epi: float = -60, fiber_space: str = 'P_1', verbose: bool = False, comm: ~mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>) Geometry[source]#

Create BiV ellipsoidal geometry

Parameters:
  • outdir (str | Path) – Directory where to save the results.

  • char_length (float, optional) – Characteristic length of mesh, by default 0.5

  • heart_as_surface (bool) – If true, create the heart as a a surface inside the torso, otherwise let the heart be a volume, by default True.

  • torso_length (float, optional) – Length of torso in the x-direction, by default 20.0

  • torso_width (float, optional) – Length of torso in the y-direction, by default 20.0

  • torso_height (float, optional) – Length of torso in the z-direction, by default 20.0

  • rotation_angle (float, optional) – Angle to rotate the torso in order to object realistic position of the heart in a torso, by default pi / 6

  • center_lv_x (float, optional) – X-coordinate for the center of the lv, by default 0.0

  • center_lv_y (float, optional) – Y-coordinate for the center of the lv, by default 0.0

  • center_lv_z (float, optional) – Z-coordinate for the center of the lv, by default 0.0

  • a_endo_lv (float, optional) – Dilation of lv endo ellipsoid in the x-direction, by default 2.5

  • b_endo_lv (float, optional) – Dilation of lv endo ellipsoid in the y-direction, by default 1.0

  • c_endo_lv (float, optional) – Dilation of lv endo ellipsoid in the z-direction, by default 1.0

  • a_epi_lv (float, optional) – Dilation of lv epi ellipsoid in the x-direction, by default 3.0

  • b_epi_lv (float, optional) – Dilation of lv epi ellipsoid in the y-direction, by default 1.5

  • c_epi_lv (float, optional) – Dilation of lv epi ellipsoid in the z-direction, by default 1.5

  • center_rv_x (float, optional) – X-coordinate for the center of the rv, by default 0.0

  • center_rv_y (float, optional) – Y-coordinate for the center of the rv, by default 0.5

  • center_rv_z (float, optional) – Z-coordinate for the center of the rv, by default 0.0

  • a_endo_rv (float, optional) – Dilation of rv endo ellipsoid in the x-direction, by default 3.0

  • b_endo_rv (float, optional) – Dilation of rv endo ellipsoid in the y-direction, by default 1.5

  • c_endo_rv (float, optional) – Dilation of rv endo ellipsoid in the z-direction, by default 1.5

  • a_epi_rv (float, optional) – Dilation of rv epi ellipsoid in the x-direction, by default 4.0

  • b_epi_rv (float, optional) – Dilation of rv epi ellipsoid in the y-direction, by default 2.5

  • c_epi_rv (float, optional) – Dilation of rv epi ellipsoid in the z-direction, by default 2.0

  • create_fibers (bool, optional) – If True create analytic fibers, by default False

  • fiber_angle_endo (float, optional) – Angle for the endocardium, by default 60

  • fiber_angle_epi (float, optional) – Angle for the epicardium, by default -60

  • fiber_space (str, optional) – Function space for fibers of the form family_degree, by default “P_1”

  • verbose (bool, optional) – If True print information from gmsh, by default False

  • comm (MPI.Comm, optional) – MPI communicator, by default MPI.COMM_WORLD

Returns:

A Geometry with the mesh, markers, markers functions and fibers.

Return type:

cardiac_geometries.geometry.Geometry

cardiac_geometries.mesh.cylinder(outdir: ~pathlib.Path | str, r_inner: float = 10.0, r_outer: float = 20.0, height: float = 40.0, char_length: float = 10.0, create_fibers: bool = False, fiber_angle_endo: float = 60, fiber_angle_epi: float = -60, fiber_space: str = 'P_1', aha: bool = True, verbose: bool = False, comm: ~mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>) Geometry[source]#

Create a cylindrical geometry

Parameters:
  • outdir (Optional[Path], optional) – Directory where to save the results.

  • r_inner (float, optional) – Radius on the endocardium layer, by default 10.0

  • r_outer (float, optional) – Radius on the epicardium layer, by default 20.0

  • height (float, optional) – Longest radius on the endocardium layer, by default 10.0

  • char_length (float, optional) – Characteristic length of mesh, by default 10.0

  • create_fibers (bool, optional) – If True create analytic fibers, by default False

  • fiber_angle_endo (float, optional) – Angle for the endocardium, by default 60

  • fiber_angle_epi (float, optional) – Angle for the epicardium, by default -60

  • fiber_space (str, optional) – Function space for fibers of the form family_degree, by default “P_1”

  • aha (bool, optional) – If True create 17-segment AHA regions

  • verbose (bool, optional) – If True print information from gmsh, by default False

  • comm (MPI.Comm, optional) – MPI communicator, by default MPI.COMM_WORLD

Returns:

A Geometry with the mesh, markers, markers functions and fibers.

Return type:

cardiac_geometries.geometry.Geometry

cardiac_geometries.mesh.cylinder_D_shaped(outdir: ~pathlib.Path | str, r_inner: float = 13.0, r_outer: float = 20.0, height: float = 40.0, inner_flat_face_distance: float = 10.0, outer_flat_face_distance: float = 17.0, char_length: float = 10.0, create_fibers: bool = False, fiber_angle_endo: float = 60, fiber_angle_epi: float = -60, fiber_space: str = 'P_1', aha: bool = True, verbose: bool = False, comm: ~mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>) Geometry[source]#

Create a D-shaped thick cylindrical shell mesh using GMSH.

Both the inner and outer surfaces have two flat faces on opposite sides.

Parameters:
  • outdir (Optional[Path], optional) – Directory where to save the results.

  • r_inner (float, optional) – Radius on the endocardium layer, by default 13.0

  • r_outer (float, optional) – Radius on the epicardium layer, by default 20.0

  • height (float, optional) – Longest radius on the endocardium layer, by default 10.0

  • inner_flat_face_distance (float) – The distance of the inner flat face from the center (along the x-axis). This value must be less than inner_radius. Default is 10.0.

  • outer_flat_face_distance (float) – The distance of the outer flat face from the center (along the x-axis). This value must be less than outer_radius. Default is 17.0.

  • char_length (float, optional) – Characteristic length of mesh, by default 10.0

  • create_fibers (bool, optional) – If True create analytic fibers, by default False

  • fiber_angle_endo (float, optional) – Angle for the endocardium, by default 60

  • fiber_angle_epi (float, optional) – Angle for the epicardium, by default -60

  • fiber_space (str, optional) – Function space for fibers of the form family_degree, by default “P_1”

  • aha (bool, optional) – If True create 17-segment AHA regions

  • verbose (bool, optional) – If True print information from gmsh, by default False

  • comm (MPI.Comm, optional) – MPI communicator, by default MPI.COMM_WORLD

Returns:

A Geometry with the mesh, markers, markers functions and fibers.

Return type:

cardiac_geometries.geometry.Geometry

cardiac_geometries.mesh.cylinder_racetrack(outdir: ~pathlib.Path | str, r_inner: float = 13.0, r_outer: float = 20.0, height: float = 40.0, inner_flat_face_distance: float = 10.0, outer_flat_face_distance: float = 17.0, char_length: float = 10.0, create_fibers: bool = False, fiber_angle_endo: float = 60, fiber_angle_epi: float = -60, fiber_space: str = 'P_1', aha: bool = True, verbose: bool = False, comm: ~mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>) Geometry[source]#

Create a racetrack-shaped thick cylindrical shell mesh using GMSH.

Both the inner and outer surfaces have two flat faces on opposite sides.

Parameters:
  • outdir (Optional[Path], optional) – Directory where to save the results.

  • r_inner (float, optional) – Radius on the endocardium layer, by default 13.0

  • r_outer (float, optional) – Radius on the epicardium layer, by default 20.0

  • height (float, optional) – Longest radius on the endocardium layer, by default 10.0

  • inner_flat_face_distance (float) – The distance of the inner flat face from the center (along the x-axis). This value must be less than inner_radius. Default is 10.0.

  • outer_flat_face_distance (float) – The distance of the outer flat face from the center (along the x-axis). This value must be less than outer_radius. Default is 17.0.

  • char_length (float, optional) – Characteristic length of mesh, by default 10.0

  • create_fibers (bool, optional) – If True create analytic fibers, by default False

  • fiber_angle_endo (float, optional) – Angle for the endocardium, by default 60

  • fiber_angle_epi (float, optional) – Angle for the epicardium, by default -60

  • fiber_space (str, optional) – Function space for fibers of the form family_degree, by default “P_1”

  • aha (bool, optional) – If True create 17-segment AHA regions

  • verbose (bool, optional) – If True print information from gmsh, by default False

  • comm (MPI.Comm, optional) – MPI communicator, by default MPI.COMM_WORLD

Returns:

A Geometry with the mesh, markers, markers functions and fibers.

Return type:

cardiac_geometries.geometry.Geometry

cardiac_geometries.mesh.lv_ellipsoid(outdir: ~pathlib.Path | str, r_short_endo: float = 7.0, r_short_epi: float = 10.0, r_long_endo: float = 17.0, r_long_epi: float = 20.0, psize_ref: float = 3, mu_apex_endo: float = -3.141592653589793, mu_base_endo: float = -1.2722641256100204, mu_apex_epi: float = -3.141592653589793, mu_base_epi: float = -1.318116071652818, create_fibers: bool = False, fiber_angle_endo: float = 60, fiber_angle_epi: float = -60, fiber_space: str = 'P_1', aha: bool = True, verbose: bool = False, comm: ~mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>) Geometry[source]#

Create an LV ellipsoidal geometry

Parameters:
  • outdir (Optional[Path], optional) – Directory where to save the results.

  • r_short_endo (float, optional) – Shortest radius on the endocardium layer, by default 7.0

  • r_short_epi (float, optional) – Shortest radius on the epicardium layer, by default 10.0

  • r_long_endo (float, optional) – Longest radius on the endocardium layer, by default 17.0

  • r_long_epi (float, optional) – Longest radius on the epicardium layer, by default 20.0

  • psize_ref (float, optional) – The reference point size (smaller values yield as finer mesh, by default 3

  • mu_apex_endo (float, optional) – Angle for the endocardial apex, by default -math.pi

  • mu_base_endo (float, optional) – Angle for the endocardial base, by default -math.acos(5 / 17)

  • mu_apex_epi (float, optional) – Angle for the epicardial apex, by default -math.pi

  • mu_base_epi (float, optional) – Angle for the epicardial apex, by default -math.acos(5 / 20)

  • create_fibers (bool, optional) – If True create analytic fibers, by default False

  • fiber_angle_endo (float, optional) – Angle for the endocardium, by default 60

  • fiber_angle_epi (float, optional) – Angle for the epicardium, by default -60

  • fiber_space (str, optional) – Function space for fibers of the form family_degree, by default “P_1”

  • aha (bool, optional) – If True create 17-segment AHA regions

  • verbose (bool, optional) – If True print information from gmsh, by default False

  • comm (MPI.Comm, optional) – MPI communicator, by default MPI.COMM_WORLD

Returns:

A Geometry with the mesh, markers, markers functions and fibers.

Return type:

cardiac_geometries.geometry.Geometry

cardiac_geometries.mesh.slab(outdir: ~pathlib.Path | str, lx: float = 20.0, ly: float = 7.0, lz: float = 3.0, dx: float = 1.0, create_fibers: bool = True, fiber_angle_endo: float = 60, fiber_angle_epi: float = -60, fiber_space: str = 'P_1', verbose: bool = False, use_dolfinx: bool = False, comm: ~mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>) Geometry[source]#

Create slab geometry

Parameters:
  • outdir (Optional[Path], optional) – Directory where to save the results.

  • lx (float, optional) – Length of slab the x-direction, by default 20.0

  • ly (float, optional) – Length of slab the x-direction, by default 7.0

  • lz (float, optional) – Length of slab the z-direction, by default 3.0

  • dx (float, optional) – Element size, by default 1.0

  • create_fibers (bool, optional) – If True create analytic fibers, by default True

  • fiber_angle_endo (float, optional) – Angle for the endocardium, by default 60

  • fiber_angle_epi (float, optional) – Angle for the epicardium, by default -60

  • fiber_space (str, optional) – Function space for fibers of the form family_degree, by default “P_1”

  • verbose (bool, optional) – If True print information from gmsh, by default False

  • use_dolfinx (bool, optional) – If True use dolfinx to create the mesh, by default False (gmsh)

  • comm (MPI.Comm, optional) – MPI communicator, by default MPI.COMM_WORLD

Returns:

A Geometry with the mesh, markers, markers functions and fibers.

Return type:

cardiac_geometries.geometry.Geometry

cardiac_geometries.mesh.slab_in_bath(outdir: ~pathlib.Path | str, lx: float = 1.0, ly: float = 0.01, lz: float = 0.5, bx: float = 0.0, by: float = 0.0, bz: float = 0.1, dx: float = 0.001, verbose: bool = False, comm: ~mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>) Geometry[source]#

Create slab geometry

Parameters:
  • outdir (Path) – Directory where to save the results.

  • lx (float, optional) – Length of slab the x-direction, by default 1.0

  • ly (float, optional) – Length of slab the x-direction, by default 0.5

  • lz (float, optional) – Length of slab the z-direction, by default 0.01

  • bx (float, optional) – Thickness of bath the x-direction, by default 0.0

  • by (float, optional) – Thickness of bath the x-direction, by default 0.0

  • bz (float, optional) – Thickness of bath the z-direction, by default 0.1

  • dx (float, optional) – Element size, by default 0.001

  • verbose (bool, optional) – If True print information from gmsh, by default False

  • comm (MPI.Comm, optional) – MPI communicator, by default MPI.COMM_WORLD

Returns:

A Geometry with the mesh, markers, markers functions and fibers.

Return type:

cardiac_geometries.geometry.Geometry

cardiac_geometries.mesh.ukb(outdir: str | ~pathlib.Path, mode: int = -1, std: float = 1.5, case: str = 'ED', char_length_max: float = 5.0, char_length_min: float = 5.0, fiber_angle_endo: float = 60, fiber_angle_epi: float = -60, create_fibers: bool = True, fiber_space: str = 'P_1', clipped: bool = False, use_burns: bool = False, burns_path: ~pathlib.Path | None = None, comm: ~mpi4py.MPI.Comm = <mpi4py.MPI.Intracomm object>) Geometry[source]#

Create a mesh from the UK-Biobank atlas using the ukb-atlas package.

Parameters:
  • outdir (str | Path) – Directory where to save the results.

  • mode (int, optional) – Mode for the UKB mesh, by default -1

  • std (float, optional) – Standard deviation for the UKB mesh, by default 1.5

  • case (str, optional) – Case for the UKB mesh (either “ED” or “ES”), by default “ED”

  • char_length_max (float, optional) – Maximum characteristic length of the mesh, by default 2.0

  • char_length_min (float, optional) – Minimum characteristic length of the mesh, by default 2.0

  • fiber_angle_endo (float, optional) – Fiber angle for the endocardium, by default 60

  • fiber_angle_epi (float, optional) – Fiber angle for the epicardium, by default -60

  • create_fibers (bool, optional) – If True create rule-based fibers, by default True

  • fiber_space (str, optional) – Function space for fibers of the form family_degree, by default “P_1”

  • clipped (bool, optional) – If True create a clipped mesh, by default False

  • use_burns (bool) – If true, use the atlas from Richard Burns to generate the surfaces. This will override the all parameter and use the burns atlas instead.

  • burns_path (Path | None) – Path to the burns atlas file. This will be a .mat file which will be loaded using scipy.io.loadmat. This needs to be specified if use_burns.

  • comm (MPI.Comm, optional) – MPI communicator, by default MPI.COMM_WORLD

Returns:

A Geometry with the mesh, markers, markers functions and fibers.

Return type:

cardiac_geometries.geometry.Geometry

class cardiac_geometries.utils.GMshGeometry(mesh, cfun, ffun, efun, vfun, markers)[source]#
cfun: MeshTags | None#

Alias for field number 1

count(value, /)#

Return number of occurrences of value.

efun: MeshTags | None#

Alias for field number 3

ffun: MeshTags | None#

Alias for field number 2

index(value, start=0, stop=9223372036854775807, /)#

Return first index of value.

Raises ValueError if the value is not present.

markers: dict[str, tuple[int, int]]#

Alias for field number 5

mesh: Mesh#

Alias for field number 0

vfun: MeshTags | None#

Alias for field number 4

class cardiac_geometries.utils.GMshModel(mesh, cell_tags, facet_tags, edge_tags, vertex_tags)[source]#
cell_tags: MeshTags#

Alias for field number 1

count(value, /)#

Return number of occurrences of value.

edge_tags: MeshTags#

Alias for field number 3

facet_tags: MeshTags#

Alias for field number 2

index(value, start=0, stop=9223372036854775807, /)#

Return first index of value.

Raises ValueError if the value is not present.

mesh: Mesh#

Alias for field number 0

vertex_tags: MeshTags#

Alias for field number 4

cardiac_geometries.utils.model_to_mesh(model, comm: ~mpi4py.MPI.Comm, rank: int, gdim: int = 3, partitioner: ~typing.Callable[[~mpi4py.MPI.Comm, int, int, ~dolfinx.cpp.graph.AdjacencyList_int32], ~dolfinx.cpp.graph.AdjacencyList_int32] | None = None, dtype=<class 'numpy.float64'>) GMshModel[source]#

Create a Mesh from a Gmsh model.

Creates a dolfinx.mesh.Mesh from the physical entities of the highest topological dimension in the Gmsh model. In parallel, the gmsh model is processed on one MPI rank, and the dolfinx.mesh.Mesh is distributed across ranks.

Parameters:
  • model – Gmsh model.

  • comm – MPI communicator to use for mesh creation.

  • rank – MPI rank that the Gmsh model is initialized on.

  • gdim – Geometrical dimension of the mesh.

  • partitioner – Function that computes the parallel distribution of cells across MPI ranks.

Returns:

A tuple containing the dolfinx.mesh.Mesh and dolfinx.mesh.MeshTags for cells, facets, edges and vertices.

Note

For performance, this function should only be called once for large problems. For re-use, it is recommended to save the mesh and corresponding tags using dolfinxio.XDMFFile after creation for efficient access.

cardiac_geometries.utils.parse_element(space_string: str, mesh: Mesh, dim: int, discontinuous: bool = False) _ElementBase[source]#

Parse a string representation of a basix element family

cardiac_geometries.utils.space_from_string(space_string: str, mesh: Mesh, dim: int, discontinuous: bool = False) 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.

  • discontinuous (bool) – If true force element to be discontinuous, by default False

class cardiac_geometries.fibers.Microstructure(f0, s0, n0)[source]#
count(value, /)#

Return number of occurrences of value.

f0: Function#

Alias for field number 0

index(value, start=0, stop=9223372036854775807, /)#

Return first index of value.

Raises ValueError if the value is not present.

n0: Function#

Alias for field number 2

s0: Function#

Alias for field number 1

cardiac_geometries.fibers.lv_ellipsoid.compute_system(t_func: Function, r_short_endo=0.025, r_short_epi=0.035, r_long_endo=0.09, r_long_epi=0.097, alpha_endo: float = -60, alpha_epi: float = 60, long_axis: int = 0, **kwargs) Microstructure[source]#

Compute the microstructure for the given time function.

Parameters:
  • t_func (dolfinx.fem.Function) – Solution of the Laplace equation

  • r_short_endo (float, optional) – Short radius at the endocardium, by default 0.025

  • r_short_epi (float, optional) – Short radius at the epicardium, by default 0.035

  • r_long_endo (float, optional) – Long radius at the endocardium, by default 0.09

  • r_long_epi (float, optional) – Long radius at the epicardium, by default 0.097

  • alpha_endo (float, optional) – Angle at the endocardium, by default -60

  • alpha_epi (float, optional) – Angle at the epicardium, by default 60

  • long_axis (int, optional) – Long axis, by default 0 (x-axis)

Returns:

The microstructure

Return type:

utils.Microstructure

cardiac_geometries.fibers.lv_ellipsoid.mu_theta(x: ndarray, y: ndarray, z: ndarray, rs: ndarray, rl: ndarray, long_axis: int = 0) tuple[ndarray, ndarray, list[int]][source]#

Get the angles mu and theta from the coordinates x, y, z given the long axis.

Parameters:
  • x (np.ndarray) – The x-coordinates

  • y (np.ndarray) – The y-coordinates

  • z (np.ndarray) – The z-coordinates

  • rs (np.ndarray) – The short radius

  • rl (np.ndarray) – The long radius

  • long_axis (int, optional) – The long axis, by default 0 (x-axis)

Returns:

The angles mu and theta and the permutation of the axes

Return type:

tuple[np.ndarray, np.ndarray, list[int]]

Raises:

ValueError – If the long axis is not 0, 1 or 2

cardiac_geometries.fibers.slab.compute_system(t_func: Function, alpha_endo: float = -60, alpha_epi: float = 60, endo_epi_axis='y') Microstructure[source]#

Compute ldrb system for slab, assuming linear angle between endo and epi

Parameters:
  • t_func (dolfin.Function) – Solution to laplace equation with 0 on endo and 1 on epi

  • alpha_endo (float, optional) – Angle on endocardium, by default -60

  • alpha_epi (float, optional) – Angle on epicardium, by default 60

Returns:

Tuple with fiber, sheet and sheet normal

Return type:

Microstructure

cardiac_geometries.fibers.slab.create_microstructure(mesh: Mesh, ffun: MeshTags, markers: dict[str, tuple[int, int]], alpha_endo: float, alpha_epi: float, function_space: str = 'P_1', outdir: str | Path | None = None, **kwargs: dict) Microstructure[source]#

Generate microstructure for slab using LDRB algorithm

Parameters:
  • mesh (dolfinx.mesh.Mesh) – A slab mesh

  • ffun (dolfinx.mesh.MeshTags) – Facet function defining the boundaries

  • markers (Dict[str, Tuple[int, int]]) – Markers with keys Y0 and Y1 representing the endo and epi planes respectively. The values should be a tuple whose first value is the value of the marker (corresponding to ffun) and the second value is the dimension

  • alpha_endo (float) – Angle on the endocardium

  • alpha_epi (float) – Angle on the epicardium

  • function_space (str) – Function space to interpolate the fibers, by default P_1

  • outdir (Optional[Union[str, Path]], optional) – Output directory to store the results, by default None. If no output directory is specified the results will not be stored, but only returned.

Returns:

Tuple with fiber, sheet and sheet normal

Return type:

Microstructure