Demo

Demo#

This script demonstrates how to generate fiber for an idealized biventricular geometry. First we import the necessary packages

import ldrb
import numpy as np
import pyvista
import dolfinx
import cardiac_geometries

We will use cardiac-geometries for generate an idealized BiV geometry and save the geometry to a folder called biv

geo = cardiac_geometries.mesh.biv_ellipsoid(outdir="biv")
2025-09-17 08:48:54 [debug    ] Convert file biv/biv_ellipsoid.msh to dolfin
Info    : Reading 'biv/biv_ellipsoid.msh'...
Info    : 43 entities
Info    : 669 nodes
Info    : 3326 elements
Info    : Done reading 'biv/biv_ellipsoid.msh'

Next we will use pyvista to plot the mesh

pyvista.start_xvfb()
vtk_mesh = dolfinx.plot.vtk_mesh(geo.mesh, geo.mesh.topology.dim)
grid = pyvista.UnstructuredGrid(*vtk_mesh)
plotter = pyvista.Plotter()
plotter.add_mesh(grid, show_edges=True)
plotter.view_xy()
if not pyvista.OFF_SCREEN:
    plotter.show()
else:
    figure = plotter.screenshot("biv_mesh.png")

It is important that the mesh is marked correctly. For a BiV geometry we need the surfaces for the endocardium (LV and RV), base and the epicardium to be marked wih a specific tag. In cardiac_geometries, the markers are saved in a dictionary

print(geo.markers)
{'BASE': [1, 2], 'ENDO_LV': [2, 2], 'ENDO_RV': [3, 2], 'EPI': [4, 2], '': [1, 3]}

Here the keys are ENDO_LV, ENDO_RV, EPI and BASE and the values are a list with two elements, the first being the marker and the second being the dimension for the marker (which is two in these cases). The ldrb algorithm accepts this format, but the default format is a dictionary with keys being lv, rv, epi and base, and the values being the values of the markers, either as an integer or as a list of integer.

Note

For single ventricular mesh, you don’t need to provide the rv key, and instead of ENDO_LV you can use ENDO

Next lets plot the facet tags with pyvista.

vtk_bmesh = dolfinx.plot.vtk_mesh(geo.mesh, geo.ffun.dim, geo.ffun.indices)
bgrid = pyvista.UnstructuredGrid(*vtk_bmesh)
bgrid.cell_data["Facet tags"] = geo.ffun.values
bgrid.set_active_scalars("Facet tags")
p = pyvista.Plotter(window_size=[800, 800])
p.add_mesh(bgrid, show_edges=True)
if not pyvista.OFF_SCREEN:
    p.show()
else:
    figure = p.screenshot("facet_tags.png")

Now let us generate fibers with 60/-60 fibers angles on the endo- and epicardium

system = ldrb.dolfinx_ldrb(
    mesh=geo.mesh,
    ffun=geo.ffun,
    markers=geo.markers,
    alpha_endo_lv=60,
    alpha_epi_lv=-60,
    beta_endo_lv=0,
    beta_epi_lv=0,
    fiber_space="P_1",
)

And let us plot the fibers with pyvista

f0 = system.f0
topology, cell_types, geometry = dolfinx.plot.vtk_mesh(f0.function_space)
values = np.zeros((geometry.shape[0], 3), dtype=np.float64)
values[:, : len(f0)] = f0.x.array.real.reshape((geometry.shape[0], len(f0)))
function_grid = pyvista.UnstructuredGrid(topology, cell_types, geometry)
function_grid["u"] = values
glyphs = function_grid.glyph(orient="u", factor=0.2)
grid = pyvista.UnstructuredGrid(*vtk_mesh)
plotter = pyvista.Plotter()
plotter.add_mesh(grid, style="wireframe", color="r")
plotter.add_mesh(glyphs)
plotter.view_xy()
if not pyvista.OFF_SCREEN:
    plotter.show()
else:
    fig_as_array = plotter.screenshot("fiber.png")