Cylindrical Geometries#

In this example we will create and visualize the three different cylindrical geometries: standard cylinder, D-shaped cylinder, and racetrack cylinder, all generated by the cardiac-geometriesx library.

from pathlib import Path
from mpi4py import MPI
import pyvista
import dolfinx
import numpy as np
import cardiac_geometries
# We disable plotting if running in a non-interactive environment, which is typical for
# automated documentation builds (e.g., Jupyter Book).
if not getattr(pyvista, "OFF_SCREEN", True):
    pyvista.OFF_SCREEN = False

1. Standard Cylindrical Shell#

The standard cylinder is a hollow cylinder defined by an inner radius, an outer radius, and a height. Analytic fibers are calculated in a circumferential-longitudinal coordinate system.

geodir_std = Path("cylinder_standard")
if not geodir_std.exists():
    cardiac_geometries.mesh.cylinder(
        outdir=geodir_std,
        char_length=10.0,
        r_inner=10.0,
        r_outer=20.0,
        height=40.0,
        create_fibers=True,
        fiber_space="P_1",
        fiber_angle_endo=60,
        fiber_angle_epi=-60,
    )
Cylindrical shell mesh generated and saved to cylinder_standard/cylinder.msh
2025-12-12 19:50:00 [debug    ] Convert file cylinder_standard/cylinder.msh to dolfin
Info    : Reading 'cylinder_standard/cylinder.msh'...
Info    : 15 entities
Info    : 169 nodes
Info    : 798 elements
Info    : Done reading 'cylinder_standard/cylinder.msh'

Mesh Visualization (Standard Cylinder)#

vtk_mesh = dolfinx.plot.vtk_mesh(geo_std.mesh, geo_std.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:
    plotter.screenshot("cylinder_standard_mesh.png")
2025-12-12 19:50:01.547 (   1.400s) [    7FD93D1A5140]vtkXOpenGLRenderWindow.:1458  WARN| bad X server connection. DISPLAY=
# The facet markers for the standard cylinder:
print(geo_std.markers)
{'INSIDE': [1, 2], 'OUTSIDE': [2, 2], 'TOP': [3, 2], 'BOTTOM': [4, 2], 'VOLUME': [5, 3]}

Facet Tags Visualization (Standard Cylinder)#

assert geo_std.ffun is not None
vtk_bmesh = dolfinx.plot.vtk_mesh(geo_std.mesh, geo_std.ffun.dim, geo_std.ffun.indices)
bgrid = pyvista.UnstructuredGrid(*vtk_bmesh)
bgrid.cell_data["Facet tags"] = geo_std.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:
    p.screenshot("cylinder_standard_facet_tags.png")

Fiber Visualization (Standard Cylinder)#

assert geo_std.f0 is not None
topology, cell_types, geometry = dolfinx.plot.vtk_mesh(geo_std.f0.function_space)
values = np.zeros((geometry.shape[0], 3), dtype=np.float64)
values[:, : len(geo_std.f0)] = geo_std.f0.x.array.real.reshape((geometry.shape[0], len(geo_std.f0)))
function_grid = pyvista.UnstructuredGrid(topology, cell_types, geometry)
function_grid["u"] = values
glyphs = function_grid.glyph(orient="u", factor=2.0)
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:
    plotter.screenshot("cylinder_standard_fiber.png")

2. D-shaped Cylinder#

The D-shaped cylinder is a thick cylindrical shell with one flat face running along its length. This geometry is created by intersecting a cylindrical shell with a cutting plane. Fiber orientations are a combination of analytic cylindrical coordinates and a slab-like orientation near the flat face, determined using the solution to the Laplace equation.

geodir_d = Path("cylinder_D_shaped")
if not geodir_d.exists():
    cardiac_geometries.mesh.cylinder_D_shaped(
        outdir=geodir_d,
        char_length=5.0,
        r_inner=13.0,
        r_outer=20.0,
        height=40.0,
        inner_flat_face_distance=10.0,
        outer_flat_face_distance=17.0,
        create_fibers=True,
        fiber_space="P_1",
    )
D-shaped cylindrical shell mesh generated and saved to cylinder_D_shaped/cylinder_D_shaped.msh
2025-12-12 19:50:02 [debug    ] Convert file cylinder_D_shaped/cylinder_D_shaped.msh to dolfin
Info    : Reading 'cylinder_D_shaped/cylinder_D_shaped.msh'...
Info    : 27 entities
Info    : 575 nodes
Info    : 2798 elements
Info    : Done reading 'cylinder_D_shaped/cylinder_D_shaped.msh'

Mesh Visualization (D-shaped Cylinder)#

vtk_mesh = dolfinx.plot.vtk_mesh(geo_d.mesh, geo_d.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:
    plotter.screenshot("cylinder_D_shaped_mesh.png")
# The facet markers for the D-shaped cylinder. Note the presence of the flat face markers (`INSIDE_FLAT` and `OUTSIDE_FLAT`):
print(geo_d.markers)
{'INSIDE_CURVED': [1, 2], 'INSIDE_FLAT': [2, 2], 'OUTSIDE_CURVED': [3, 2], 'TOP': [4, 2], 'OUTSIDE_FLAT': [5, 2], 'BOTTOM': [6, 2], 'VOLUME': [7, 3]}

Facet Tags Visualization (D-shaped Cylinder)#

assert geo_d.ffun is not None
vtk_bmesh = dolfinx.plot.vtk_mesh(geo_d.mesh, geo_d.ffun.dim, geo_d.ffun.indices)
bgrid = pyvista.UnstructuredGrid(*vtk_bmesh)
bgrid.cell_data["Facet tags"] = geo_d.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:
    p.screenshot("cylinder_D_shaped_facet_tags.png")

Fiber Visualization (D-shaped Cylinder)#

assert geo_d.f0 is not None
topology, cell_types, geometry = dolfinx.plot.vtk_mesh(geo_d.f0.function_space)
values = np.zeros((geometry.shape[0], 3), dtype=np.float64)
values[:, : len(geo_d.f0)] = geo_d.f0.x.array.real.reshape((geometry.shape[0], len(geo_d.f0)))
function_grid = pyvista.UnstructuredGrid(topology, cell_types, geometry)
function_grid["u"] = values
glyphs = function_grid.glyph(orient="u", factor=2.0)
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:
    plotter.screenshot("cylinder_D_shaped_fiber.png")

3. Racetrack Cylinder#

The racetrack cylinder is a variant of the D-shaped cylinder with two parallel flat faces. This geometry is created by intersecting the cylindrical shell with two parallel cutting planes. The fiber generation method combines the analytic cylindrical model with two slab-like regions at the flat faces.

geodir_race = Path("cylinder_racetrack")
if not geodir_race.exists():
    cardiac_geometries.mesh.cylinder_racetrack(
        outdir=geodir_race,
        char_length=5.0,
        r_inner=13.0,
        r_outer=20.0,
        height=40.0,
        inner_flat_face_distance=10.0,
        outer_flat_face_distance=17.0,
        create_fibers=True,
        fiber_space="P_1",
    )
D-shaped cylindrical shell mesh generated and saved to cylinder_racetrack/cylinder_racetrack.msh
2025-12-12 19:50:03 [debug    ] Convert file cylinder_racetrack/cylinder_racetrack.msh to dolfin
Info    : Reading 'cylinder_racetrack/cylinder_racetrack.msh'...
Info    : 51 entities
Info    : 629 nodes
Info    : 3071 elements
Info    : Done reading 'cylinder_racetrack/cylinder_racetrack.msh'

Mesh Visualization (Racetrack Cylinder)#

vtk_mesh = dolfinx.plot.vtk_mesh(geo_race.mesh, geo_race.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:
    plotter.screenshot("cylinder_racetrack_mesh.png")
# The facet markers for the racetrack cylinder. Note the markers for the two sets of inner/outer flat faces (`INSIDE_FLAT1`, `OUTSIDE_FLAT1`, `INSIDE_FLAT2`, `OUTSIDE_FLAT2`):
print(geo_race.markers)
{'INSIDE_CURVED1': [1, 2], 'INSIDE_FLAT1': [2, 2], 'INSIDE_FLAT2': [3, 2], 'INSIDE_CURVED2': [4, 2], 'OUTSIDE_CURVED1': [5, 2], 'OUTSIDE_FLAT1': [6, 2], 'BOTTOM': [7, 2], 'OUTSIDE_FLAT2': [8, 2], 'TOP': [9, 2], 'OUTSIDE_CURVED2': [10, 2], 'VOLUME': [11, 3]}

Facet Tags Visualization (Racetrack Cylinder)#

assert geo_race.ffun is not None
vtk_bmesh = dolfinx.plot.vtk_mesh(geo_race.mesh, geo_race.ffun.dim, geo_race.ffun.indices)
bgrid = pyvista.UnstructuredGrid(*vtk_bmesh)
bgrid.cell_data["Facet tags"] = geo_race.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:
    p.screenshot("cylinder_racetrack_facet_tags.png")

Fiber Visualization (Racetrack Cylinder)#

assert geo_race.f0 is not None
topology, cell_types, geometry = dolfinx.plot.vtk_mesh(geo_race.f0.function_space)
values = np.zeros((geometry.shape[0], 3), dtype=np.float64)
values[:, : len(geo_race.f0)] = geo_race.f0.x.array.real.reshape((geometry.shape[0], len(geo_race.f0)))
function_grid = pyvista.UnstructuredGrid(topology, cell_types, geometry)
function_grid["u"] = values
glyphs = function_grid.glyph(orient="u", factor=2.0)
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:
    plotter.screenshot("cylinder_racetrack_fiber.png")