svVascularize

End-to-End Pipeline for 3D CFD Simulations

Generating a watertight, simulation-ready finite element mesh for computational fluid dynamics is a challenging multi-step problem. For typical workflows modeling bloodflow from clinical imaging data, this process involves pathline generation, segmentation, solid surface creation, volume mesh generation, mesh optimization/smoothing, and boundary identification. The pathline and segmentation steps are dependent upon available image data, thus these steps are very different for synthetically generated vascular structures which inheriently lack imaging data. The remaining steps are somewhat the same except that the workflow employed with svVascularize is automated since the number of vessels may quickly exceed practical manual efforts for model and mesh creation. These automated steps are all part of the Simulation class which can be used to wrap around synthetic vascular objects in order to build simulation files. We will demonstrate the basic API below for a simple vascular tree.

 import pyvista as pv
 from svv.domain.domain import Domain
 from svv.tree.tree import Tree
 from svv.simulation.simulation import Simulation

 cube = Domain(pv.Cube())
 cube.create()
 cube.solve()
 cube.build()

 t = Tree()
 t.set_domain(cube)
 t.set_root()
 t.n_add(2)

Now that we have a simple tree we can start building simulations from this object.

 # Create Simulation Container for Synthetic Tree Object
 sim = Simulation(t)

 # Build All Surface/Volume Meshes for 3D CFD Simulation
 sim.build_meshes()

This method will produce unioned surface .vtp meshes and their cooresponding volume meshes .vtu representing the synthetic vascular object. Additionally, the procedure will attempt to append boundary layer meshes along the walls vessels, typical in many hemodynamic studies to account for the steep velocity gradients near vessel walls. If successful, the surfaces meshes will be stored within the simulation container list simulation.fluid_domain_surface_meshes, volume meshes will be stored within the container simulation.fluid_domain_volume_meshes, and boundary layers will be stored within the container simulation.fluid_domain_boundary_layers. Storing the meshes as python objects allows for more custom modification procedures that users may want to implement prior to creating the simulation files.

Mesh building: This step is potentially memory and compute intensive depending on the size of the synthetic vascular object being discretized. Under the hood, svVascularize is using a combination of tetgen and mmg meshing utilities via python subprocess calls. These, as well as surface mesh unioning operations, are non-trivial and may take minutes-hours to complete. Users can supply custom meshing parameters for their applications to override default values this may result in more efficient meshing, but it is recommended that default parameters are attemped first.

Once meshes exist, identify wall and outlet patches that become boundary conditions. The Simulation.extract_faces() routine wraps svv.simulation.utils.extract_faces to label surfaces on both the solid model and the tetrahedral mesh.

 # Detect walls, outlets, and shared boundaries
sim.extract_faces(crease_angle=60.0)

# Inspect names stored in the GeneralMesh container
print(sim.fluid_domain_meshes[0].faces.keys())

Finally, assemble the svFSI configuration and write everything to disk. The default output directory is created if it does not exist.

 # Build and persist 3-D CFD input files
sim.construct_3d_fluid_simulation()
sim.write_3d_fluid_simulation()

# Result: meshes + XML under sim.file_path (e.g. ./simulations/example/)
Directory layout: The export creates a mesh/ folder containing the volume grid (.vtu or .vtp), a mesh-surfaces/ folder for each named face, and a fluid_simulation_0-0.xml svFSI input file.

For rapid prototyping, couple the same tree to the reduced-order solvers shipped with SimVascular.

 from svv.simulation.fluid.rom.zero_d.zerod_tree import export_0d_simulation

export_0d_simulation(t, outdir=sim.file_path + "/rom_0d")

# Optionally stage 1-D parameters for multi-scale workflows
sim.construct_1d_fluid_simulation(time_step_size=1e-3, number_time_steps=2000)
Note: Writing of 1-D files is not yet automated. The parameter bundle returned by construct_1d_fluid_simulation() can be passed to the SimVascular 1-D solver tools.