API Reference

Contents

API Reference#

Pre-processing scripts#

vampy.automatedPreprocessing.automated_preprocessing.main_meshing()[source]#
vampy.automatedPreprocessing.automated_preprocessing.read_command_line(input_path=None)[source]#

Read arguments from commandline and return all values in a dictionary. If input_path is not None, then do not parse command line, but only return default values.

Args:

input_path (str): Input file path, positional argument with default None.

vampy.automatedPreprocessing.automated_preprocessing.run_pre_processing(input_model, verbose_print, smoothing_method, smoothing_factor, smoothing_iterations, meshing_method, refine_region, is_atrium, add_flow_extensions, visualize, config_path, coarsening_factor, inlet_flow_extension_length, outlet_flow_extension_length, edge_length, region_points, compress_mesh, add_boundary_layer, scale_factor, resampling_step, flow_rate_factor, moving_mesh, clamp_boundaries, max_geodesic_distance)[source]#

Automatically generate mesh of surface model in .vtu and .xml format, including prescribed flow rates at inlet and outlet based on flow network model.

Runs simulation of meshed case on a remote ssh server if server configuration is provided.

Parameters
  • max_geodesic_distance (float) – Max distance when performing geodesic refinement (in mm)

  • input_model (str) – Name of case

  • verbose_print (bool) – Toggles verbose mode

  • smoothing_method (str) – Method for surface smoothing

  • smoothing_factor (float) – Smoothing factor of Voronoi smoothing

  • smoothing_iterations (int) – Number of smoothing iterations for Taubin and Laplace smoothing

  • meshing_method (str) – Determines what the density of the volumetric mesh depends upon

  • refine_region (bool) – Refines selected region of input if True

  • is_atrium (bool) – Determines whether this is an atrium case

  • add_flow_extensions (bool) – Adds flow extensions to mesh if True

  • visualize (bool) – Visualize resulting surface model with flow rates

  • config_path (str) – Path to configuration file for remote simulation

  • coarsening_factor (float) – Refine or coarsen the standard mesh size with given factor

  • region_points (list) – User defined points to define which region to refine

  • edge_length (float) – Edge length used for meshing with constant element size

  • inlet_flow_extension_length (float) – Factor defining length of flow extensions at the inlet(s)

  • outlet_flow_extension_length (float) – Factor defining length of flow extensions at the outlet(s)

  • compress_mesh (bool) – Compresses finalized mesh if True

  • add_boundary_layer (bool) – Adds boundary layers to walls if True

  • scale_factor (float) – Scale input model by this factor

  • resampling_step (float) – Float value determining the resampling step for centerline computations, in [m]

  • flow_rate_factor (float) – Flow rate factor

  • moving_mesh (bool) – Computes projected movement for displaced surfaces located in [filename_model]_moved folder

  • clamp_boundaries (bool) – Clamps inlet(s) and outlet(s) if true

vampy.automatedPreprocessing.preprocessing_common.add_flow_extension(surface, centerlines, is_inlet, extension_length=2.0, extension_mode='boundarynormal')[source]#

Adds flow extensions to either all inlets or all outlets with specified extension length.

Parameters
  • surface (vtkPolyData) – Surface model to extend

  • centerlines (vtkPolyData) – Centerlines in model.

  • is_inlet (bool) – Determines if inlet or outlet is to be extended.

  • extension_length (float) – Determines length of flow extensions. Factor is multiplied with MISR at relevant boundary

  • extension_mode (str) – Determines how extensions are place, either normal to boundary or following centerline direction

Returns

Extended surface model

Return type

surface_extended (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.check_if_closed_surface(surface)[source]#

Checks if the given surface is capped (i.e., has no feature edges).

Parameters

surface (vtkPolyData) – The surface to check for capping.

Returns

True if the surface is capped, False otherwise.

Return type

bool

vampy.automatedPreprocessing.preprocessing_common.compute_centers_for_meshing(surface, is_atrium, case_path=None, test_capped=False)[source]#

Compute the center of all the openings in the surface. The inlet is chosen based on the largest area for arteries (or aneurysm). However, for atrium, the outlet is chosen based on the largest area (new).

Parameters
  • test_capped (bool) – Check if surface is capped.

  • surface (vtkPolyData) – centers of the openings.

  • case_path (str) – path to case directory.

  • is_atrium (bool) – Check if it is an atrium model.

Returns

Inlet center. outlet_centers (list): A flattened list with all the outlet centers.

Return type

inlet_center (list)

vampy.automatedPreprocessing.preprocessing_common.compute_distance_to_sphere(surface, center_sphere, radius_sphere=0.0, distance_offset=0.0, distance_scale=0.01, min_distance=0.2, max_distance=0.3)[source]#

Computes cell specific target edge length (distances) based on input criterion.

Parameters
  • surface (vtkPolyData) – Input surface model

  • center_sphere (list) – Point representing region to refine

  • radius_sphere (list) – Maximum inscribed sphere radius in region of refinement

  • distance_offset (float) – Offsets the relevant region by a offset

  • distance_scale (float) – Determines how the distance is scaled based on the distance from the relevant region

  • min_distance (float) – Minimum distance away from the relevant region before scaling starts

  • max_distance (float) – Maximum distance away from the relevant region before scaling stops

Returns

Modified surface model with distances

Return type

surface (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.compute_flow_rate(is_atrium, inlet, parameters, flow_rate_factor)[source]#

Computes mean flow rate used as a boundary condition for pulsatile flow condition

Parameters
  • is_atrium (bool) – Determines if model is atrium or artery

  • inlet (list) – List of points representing midpoint of boundaries

  • parameters (dict) – Dictionary containing model parameters

  • flow_rate_factor (float) – Factor to adjust flow rate calculation

Returns

Mean inflow rate

Return type

mean_inflow_rate (float)

vampy.automatedPreprocessing.preprocessing_common.dist_sphere_constant(surface, centerlines, region_center, misr_max, save_path, edge_length)[source]#

Determines the target edge length for each cell on the surface, including potential refinement or coarsening of certain user specified areas. Level of refinement/coarseness is determined based on user selected region, otherwise a constant target edge length is selected.

Parameters
  • surface (vtkPolyData) – Input surface model

  • centerlines (vtkPolyData) – Centerlines of input model

  • region_center (list) – Point representing region to refine

  • misr_max (list) – Maximum inscribed sphere radius in region of refinement

  • save_path (str) – Location to store processed surface

  • edge_length (float) – Target edge length

Returns

Processed surface model with info on cell specific target edge length

Return type

surface (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.dist_sphere_curvature(surface, centerlines, region_center, misr_max, save_path, factor)[source]#

Determines the target edge length for each cell on the surface, including potential refinement or coarsening of certain user specified areas. Level of refinement/coarseness is determined based on surface curvature.

Parameters
  • surface (vtkPolyData) – Input surface model

  • centerlines (vtkPolyData) – Centerlines of input model

  • region_center (list) – Point representing region to refine

  • misr_max (list) – Maximum inscribed sphere radius in region of refinement

  • save_path (str) – Location to store processed surface

  • factor (float) – Coarsening factor, determining the level of refinement (<1) or coarsening (>1)

Returns

Processed surface model with info on cell specific target edge length

Return type

surface (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.dist_sphere_diam(surface, centerlines, region_center, misr_max, save_path, factor)[source]#

Determines the target edge length for each cell on the surface, including potential refinement or coarsening of certain user specified areas. Level of refinement/coarseness is determined based on the distance to the centerline.

Parameters
  • surface (vtkPolyData) – Input surface model

  • centerlines (vtkPolyData) – Centerlines of input model

  • region_center (list) – Point representing region to refine

  • misr_max (list) – Maximum inscribed sphere radius in region of refinement

  • save_path (str) – Location to store processed surface

  • factor (float) – Coarsening factor, determining the level of refinement (<1) or coarsening (>1)

Returns

Processed surface model with info on cell specific target edge length

Return type

surface (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.dist_sphere_geodesic(surface, region_center, max_distance, save_path, edge_length)[source]#

Determines the target edge length for each cell on the surface, including potential refinement or coarsening of certain user specified areas. Level of refinement/coarseness is determined based on user selected region and the geodesic distance from it. :param surface: Input surface model :type surface: vtkPolyData :param region_center: Point representing region to refine :type region_center: list :param save_path: Location to store processed surface :type save_path: str :param edge_length: Target edge length :type edge_length: float :param max_distance: Max distance of geodesic measure :type max_distance: float

Returns

Processed surface model with info on cell specific target edge length

Return type

surface (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.find_boundaries(model_path, mean_inflow_rate, network, mesh, verbose_print, is_atrium)[source]#

Finds inlet and outlet boundary IDs after complete meshing, including mean flow ratio and area ratios between outlets or inlets (determined by type of model)

Parameters
  • model_path (str) – Path to model files

  • mean_inflow_rate (float) – Flow rate

  • network (Network) – Flow splitting network based on network boundary condition

  • mesh (vtkUnstructuredGrid) – Volumetric mesh

  • verbose_print (bool) – Prints additional info if True

  • is_atrium (bool) – Determines if model represents atrium or artery

vampy.automatedPreprocessing.preprocessing_common.generate_mesh(surface, add_boundary_layer)[source]#

Generates a mesh suitable for CFD from a input surface model.

Parameters

surface (vtkPolyData) – Surface model to be meshed.

Returns

Output mesh remeshedsurface (vtkPolyData): Remeshed version of the input model

Return type

mesh (vtkUnstructuredGrid)

vampy.automatedPreprocessing.preprocessing_common.get_centers_for_meshing(surface, is_atrium, dir_path, use_flow_extensions=False)[source]#

Get the centers of the inlet and outlets.

Parameters
  • surface (vtkPolyData) – An open surface.

  • is_atrium (bool) – True if model represents atrium

  • dir_path (str) – Path to the case file.

  • use_flow_extensions (bool) – Turn on/off flow extension.

Returns

A flatt list with the point of the inlet outlet (list): A flatt list with the points of all the outlets.

Return type

inlet (list)

vampy.automatedPreprocessing.preprocessing_common.get_furtest_surface_point(inlet, surface)[source]#

Calculates the furthest point on a given surface from a specified inlet point. It iterates over each point on the surface, calculating its distance from the inlet and updates the furthest distance accordingly.

Parameters
  • inlet (list) – A list with the inlet point.

  • surface (vtkPolyData) – The surface to check for the furthest point.

Returns

The x, y, z coordinates of the furthest point on the surface from the inlet.

Return type

outlets (list)

vampy.automatedPreprocessing.preprocessing_common.get_regions_to_refine(surface, provided_points, dir_path)[source]#

Determines which regions to refine.

Parameters
  • surface (vtkPoylData) – Surface model

  • provided_points (ndarray) – Points defining regions to refine

  • dir_path (str) – Path to save directory

Returns

List of points representing regions to refine

Return type

region_points (list)

vampy.automatedPreprocessing.preprocessing_common.make_voronoi_diagram(surface, file_path)[source]#

Compute the voronoi diagram of surface model.

Parameters
  • surface (polydata) – Capped surface model to create a Voronoi diagram of.

  • file_path (str) – Absolute path to surface model path.

Returns

Voronoi diagram of surface.

Return type

voronoi (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.mesh_alternative(surface)[source]#

Subdivides and smoothes the input surface model, preparing it for volumetric meshing.

Parameters

surface (vtkPolyData) – Input surface model to be meshed alternatively

Returns

Smoothed model

Return type

surface (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.provide_region_points(surface, provided_points, dir_path=None)[source]#

Get relevant region points from user selected points on a input surface.

Parameters
  • provided_points (ndarray) – Point(s) representing area to refine

  • surface (vtkPolyData) – Surface model.

  • dir_path (str) – Location of info.json file

Returns

List of relevant outlet IDs

Return type

points (list)

vampy.automatedPreprocessing.preprocessing_common.remesh_surface(surface, edge_length, element_size_mode='edgelength', exclude=None)[source]#

Remeshes a given surface based on the specified parameters.

Parameters
  • surface (vtkPolyData) – The input surface to be remeshed.

  • edge_length (float) – The target edge length for remeshing.

  • element_size_mode (str, optional) – Determines the method for sizing elements during remeshing.

  • exclude (list, optional) – A list of entity IDs to be excluded during remeshing. Defaults to None.

Returns

The remeshed surface.

Return type

remeshed_surface (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.scale_surface(surface, scale_factor)[source]#

Scale the input surface by a factor scale_factor.

Parameters
  • surface (vtkPolyData) – Input surface to be scaled

  • scale_factor (float) – Scaling factor

Returns

Scaled input surface

Return type

scaled_surface (vtkPolyData)

vampy.automatedPreprocessing.preprocessing_common.setup_model_network(centerlines, file_name_probe_points, region_center, verbose_print, is_atrium)[source]#

Sets up network used for network boundary condition model.

Parameters
  • centerlines (vtkPolyData) – Centerlines representing meshed model

  • file_name_probe_points (str) – Save path of probe points

  • region_center (list) – List of points representing region of refinement

  • verbose_print (bool) – Prints additional info if True

  • is_atrium (bool) – Determines if model is atrium or artery

Returns

Network model probe_points (ndarray): Probe points where velocity and pressure is to be sampled

Return type

network (Network)

vampy.automatedPreprocessing.preprocessing_common.write_mesh(compress_mesh, file_name_surface_name, file_name_vtu_mesh, file_name_xml_mesh, mesh, remeshed_surface)[source]#

Writes the mesh to DOLFIN format, and compresses to .gz format

Parameters
  • compress_mesh (bool) – Compressed mesh to zipped format

  • file_name_surface_name (str) – Path to remeshed surface model

  • file_name_vtu_mesh (str) – Path to VTK mesh

  • file_name_xml_mesh (str) – Path to XML mesh

  • mesh (vtuUnstructuredGrid) – Meshed surface model

  • remeshed_surface (vtkPolyData) – Remeshed surface model

vampy.automatedPreprocessing.moving_common.easing_cos(dist_b, dist_bb)[source]#

A simple easing function between 0 and 1 for creating a sinusoidal profile :param dist_b: Distance to boundary for current point :type dist_b: ndarray :param dist_bb: Distance between boundaries :type dist_bb: ndarray

Returns

Sinusoidal profile based on distance between point and boundary

Return type

distances (ndarray)

vampy.automatedPreprocessing.moving_common.get_point_map(remeshed, remeshed_extended, profile='linear')[source]#

Create a map between original surface model and model with flow extensions; remeshed_extended, for clamping the inlet(s) and outlet. A linear or sinusodial profile for movement between endpoints and moving domain is prescribed. :param remeshed: Input surface :type remeshed: vtkPolyData :param remeshed_extended: Input surface with flow extensions :type remeshed_extended: vtkPolyData :param profile: ‘linear’ or ‘sine’ for respective profiles :type profile: str

Returns

Array of linear or sinusodial profile between boundary and moving mesh point_map (ndarray): Mapping between surface IDs and points along extension

Return type

distances (ndarray)

vampy.automatedPreprocessing.moving_common.move_surface_model(surface, original, remeshed, remeshed_extended, distance, point_map, file_path, i, points, clamp_boundaries)[source]#

Computes and applies the displacement between the original and the given surface, then projects this displacement onto a remeshed surface. The displaced points of the remeshed surface are then saved to a specified file path.

Parameters
  • surface (vtkPolyData) – The source surface used to compute the displacement.

  • original (vtkPolyData) – The reference surface for computing the displacement.

  • remeshed (vtkPolyData) – A remeshed version of the original surface.

  • remeshed_extended (vtkPolyData) – An extended version of the remeshed surface.

  • distance (float) – Scalar factor for the displacement.

  • point_map (ndarray) – Mapping of the points from remeshed to remeshed_extended.

  • file_path (str) – Path to save the displaced points of the remeshed surface.

  • i (int) – Index to specify where in the ‘points’ array the data should be saved.

  • points (ndarray) – An array where the displaced points will be saved.

  • clamp_boundaries (bool) – If True, the displacement will be clamped by the distance in the extensions.

Returns

The function writes results to a file and modifies the points array in place.

Return type

None

vampy.automatedPreprocessing.moving_common.project_displacement(clamp_boundaries, distance, folder_extended_surfaces, folder_moved_surfaces, point_map, surface, surface_extended, remeshed, scale_factor)[source]#

Projects the displacement of moved surfaces located in a specified folder onto an extended surface. The projected surfaces are then saved to a designated output directory. Only processes files in “folder_moved_surfaces” that end with “.vtp” or “.stl”.

Parameters
  • clamp_boundaries (bool) – If True, the displacement will be clamped by the distance in the extensions.

  • distance (float) – Scalar factor for the displacement.

  • folder_extended_surfaces (str) – Path to the directory where extended surfaces will be saved.

  • folder_moved_surfaces (str) – Path to the directory containing the moved surfaces to be projected.

  • point_map (ndarray) – Mapping of the points from remeshed to surface_extended.

  • surface (vtkPolyData) – The reference surface for computing the displacement.

  • surface_extended (vtkPolyData) – An extended version of the reference surface.

  • remeshed (vtkPolyData) – A remeshed version of the original surface.

  • scale_factor (float) – Factor to scale the surfaces.

Returns

A 3D numpy array where the displaced points will be saved. The shape is (num_points, 3, num_surfaces).

Return type

array

vampy.automatedPreprocessing.moving_common.save_displacement(file_name_displacement_points, points, number_of_points=100)[source]#

Resample displacement points and write them to numpy data array :param file_name_displacement_points: Path to numpy point path :type file_name_displacement_points: str :param points: Numpy array consisting of displacement surface points :type points: ndarray :param number_of_points: Number of points to interpolate displacement points over :type number_of_points: int

class vampy.automatedPreprocessing.DisplayData.DisplayModel[source]#

Bases: object

DisplayProbesAndModel(centerline, fileNameCenterline, listProbePoints, model=None)[source]#

Displays a model and the corresponding probe points along the centerline.

polyDataToActor(polyData, opacity=1.0)[source]#

Wrap the provided vtkPolyData object in a mapper and an actor, returning the actor.

renderWindow(renderer, titleWindow)[source]#
setLight(renderer)[source]#
class vampy.automatedPreprocessing.DisplayData.VtkPointCloud(maxNumPoints=1000000.0)[source]#

Bases: object

addPoint(point)[source]#
clearPoints()[source]#
vampy.automatedPreprocessing.ImportData.ComputeBranchLength(centerline, branchId)[source]#

Return the length for the branch of index branchId.

vampy.automatedPreprocessing.ImportData.ComputeBranchRadius(centerline, branchId)[source]#

Return the radius for a branch with index branchId.

vampy.automatedPreprocessing.ImportData.ComputeConnectivity(network, tolerance, verboseprint)[source]#

Compute the branches connectivity in the network.

It is assumed that the first element is the inlet. In case of several inlets, modification to this routine should be coded. This routine search for each branch, its connected branches by comparing the coordinates of the distal extremity of the treated branch with the proximal extremity of the other branches in the network. As sometimes the coordinates do not perfectly match, a tolerance is applied.

TODO: use simply vtkPolyDataConnectivityFilter to do the job…?

vampy.automatedPreprocessing.ImportData.ComputeGeometricTolerance(centerline)[source]#

Return the min and max length for branches.

This routine compute the delta x minimum and maximum in the network. The minimum length is used for the connectivity. The maximum length is used for merging the blanked segments.

vampy.automatedPreprocessing.ImportData.ComputeGroupLength(centerline, branchGroupId)[source]#

Return the mean length for branches of branchGroupId.

vampy.automatedPreprocessing.ImportData.ComputeGroupRadius(centerline, branchGroupId)[source]#

Return the mean radius of a group.

The mean radius is computed using the hydraulic resistance for the branches of branchGroupId. See Chnafa et al. International Journal for Numerical Methods in Biomedical Engineering, 2016.

vampy.automatedPreprocessing.ImportData.ComputeInletAverageCrossSectionArea(centerline, desiredLength, verboseprint)[source]#

Compute the inlet radius as an averaged radius.

Computes an average radius over a certain length of the ICA. The mean radius is computed with the cross section area. If the number of points over the segment is lower than two, then the section area correspond to the inlet point. Note that due to a bug in VTK the very first point seems wrong, that is why I am taking the second point of the centerline.

vampy.automatedPreprocessing.ImportData.ComputeInletAverageRadius(centerline, desiredLength, verboseprint)[source]#

Compute the inlet radius as an averaged radius.

Computes an average radius over a certain length of the ICA. The mean radius is computed with the maximum inscribed sphere radius in the vessel. Thus, the actual radius might be underestimated if the vessel has an elliptical shape.

class vampy.automatedPreprocessing.ImportData.Element(Id)[source]#

Bases: object

GetAlpha()[source]#
GetBehindSegment()[source]#
GetBeta()[source]#
GetFrontSegment()[source]#
GetGamma()[source]#
GetId()[source]#
GetInPointsx0()[source]#
GetInPointsx0Id()[source]#
GetInletRadius()[source]#
GetLength()[source]#
GetMeanArea()[source]#
GetMeanRadius()[source]#
GetOutPointsx1()[source]#
GetOutPointsx1Id()[source]#
GetVtkCellIdList()[source]#
GetVtkGroupIdList()[source]#
IsAnInlet()[source]#
IsAnOutlet()[source]#
IsBlanked()[source]#

Returns True if the element is part of a branch division.

SetAlpha(alpha)[source]#
SetBehindSegment(id)[source]#
SetBeta(beta)[source]#
SetBlanking(blanking)[source]#
SetFrontSegment(id)[source]#
SetGamma(gamma)[source]#
SetIfInlet(inlet)[source]#
SetIfOutlet(outlet)[source]#
SetInOutPointsCoordinates(x0, x1)[source]#
SetInOutPointsIds(x0Id, x1Id)[source]#
SetInletRadius(radius)[source]#
SetLength(length)[source]#
SetMeanRadius(radius)[source]#
SetOutletRadius(radius)[source]#
SetVtkCellIdList(vtkCellIdList)[source]#
SetVtkGroupIdList(VtkGroupIdList)[source]#
vampy.automatedPreprocessing.ImportData.GetBlankedGroupsIdList(centerline)[source]#
vampy.automatedPreprocessing.ImportData.GetIndexCenterlineForADefinedLength(centerline, branchId, desiredLength, verboseprint)[source]#

Get the index of the point such as the desired distance between the index and the beginning of the branch is reached.

vampy.automatedPreprocessing.ImportData.GetListProbePoints(centerline, network, verboseprint)[source]#

Get points on a centerline spaced by the inlet diameter.

vampy.automatedPreprocessing.ImportData.GetListsUniqueBlankedBranches(blankedGroupsIdList, redundantBlankedBranchesIdList)[source]#
vampy.automatedPreprocessing.ImportData.GetMaxGroupId(centerline)[source]#
vampy.automatedPreprocessing.ImportData.GetRedundantBlankedIdList(centerline, blankedGroupsIdList)[source]#

Get the redundant blanked segments.

Nominally, the x0 and x1 coordinates should be the same. Thus, the computed distance should be exactly 0.0. Sometimes branches are just different by a few and SHOULD be merged. Visual inspection of the centerline with Paraview when the Warning message pop up is advised. Sometimes, the problem is more concerning. It is problematic that VMTK split, sometimes, the branches with a segment with a blanked segment, a lil’ something and another blanked segment. In case of wierd branch splittings, try to lower the tol. It is cheating since branches that should not be merge will be merged anyway.

vampy.automatedPreprocessing.ImportData.IsArrayDefined(centerline, arrayName)[source]#
class vampy.automatedPreprocessing.ImportData.Network[source]#

Bases: object

AddElement(x)[source]#
GetNetworkInletRadius()[source]#
GetNumberOfBifBranches()[source]#
GetNumberOfElements()[source]#
GetNumberOfInlet()[source]#
GetNumberOfOutlet()[source]#
SetNetworkInletRadius(radius)[source]#
vampy.automatedPreprocessing.ImportData.SetNetworkStructure(centerline, network, verboseprint, isConnectivityNeeded=True, isRadiusInletNeeded=True)[source]#

Fills a network structure with a vtkPolyData object.

Each element has an unique index. The groups length and radius are averaged out. When the element is not not part of the bifurcation reference system, it is possible that multiple end points exist. They are all saved. For the complet formalism of the bifurcation reference system see Antiga, L., & Steinman, D. A. (2004). Robust and objective decomposition and mapping of bifurcating vessels. Medical Imaging, IEEE Transactions on, 23(6), 704-713.

vampy.automatedPreprocessing.ImportData.SetRadiusX0(centerline, network, verboseprint)[source]#
vampy.automatedPreprocessing.ImportData.loadFile(fileName)[source]#

Load the given file, and return a vtkPolyData object for it.

class vampy.automatedPreprocessing.NetworkBoundaryConditions.FlowSplitting[source]#

Bases: object

This class computes the flow splitting for a network according to [ADD REF AJRN] and [VALIDATION PAPER]

CheckTotalFlowRate(network, verboseprint)[source]#

Check if the sum of the outflows is 100%.

ComputeAlphas(network, verboseprint)[source]#

Compute the segments alpha coefficient.

Compute for each blanked segment of a network its attributed alpha coefficient. For a network division with N daughter segments (N can be superior to two), the alpha coefficient of each segment i is computed as: alpha_i = S / sumSurfaces = S_i / sum_1^N(segment area_j) That means that alpha is the percentage of flow from the datum (the supplier) segment going though the i-th segment of the considered division. Segments that are not blanked are left with an alpha coefficient of 1.0 since there is no flow division.

ComputeBetas(network, verboseprint)[source]#

Compute the outlets beta coefficient.

Compute for each outlet segment of a network its attributed beta coefficient. This coefficient correspond to the percentage of flow going out of this segment from the inlet flow. It is the actual flow division at the outlets. Segments that are not outlets are left with a beta coefficient of 0.0. The algorithm compute the beta coefficients as the multiplications of the alpha coefficient from the considered outlet to the root of the network where the inlet is applied.

ComputeGammas(network, verboseprint)[source]#

Compute the outlets gamma coefficient.

Compute for each outlet segment of a network its attributed beta coefficient. This coefficient correspond to the percentage of flow going out of this segment from the inlet flow. The Gamma coefficient of each outlet i is computed as: gamma_i = S_i / sumOutletAreas

vampy.automatedPreprocessing.repair_tools.clean_surface(surface)[source]#
vampy.automatedPreprocessing.repair_tools.find_and_delete_nan_triangles(surface)[source]#
vampy.automatedPreprocessing.repair_tools.is_nan(num)[source]#
vampy.automatedPreprocessing.repair_tools.print_surface_info(surface)[source]#
vampy.automatedPreprocessing.simulate.exists(sftp, path)[source]#

os.path.exists for paramiko’s SCP object

vampy.automatedPreprocessing.simulate.run_simulation(config_path, local_dir, case_name)[source]#

Run simulation of case on a remote ssh server with given input configuration.

Parameters
  • config_path (str) – Path to configuration file

  • local_dir (str) – Path to case folder

  • case_name (str) – Case name

vampy.automatedPreprocessing.visualize.visualize_model(network_elements, probe_points, output_surface, mean_inflow_rate)[source]#

Visualize surface model / mesh with distributed flow rates at each inlet/outlet in an interactive VTK window.

Parameters
  • network_elements (list) – List of inlet/outlet elements within a surface model network

  • probe_points (ndarray) – Array of probe points

  • output_surface (vtkPolyData) – Surface model to visualize

  • mean_inflow_rate (float) – Mean inflow rate of input model

class vampy.automatedPreprocessing.vmtk_pointselector.VtkText(guiText='')[source]#

Bases: object

class vampy.automatedPreprocessing.vmtk_pointselector.vmtkPickPointSeedSelector[source]#

Bases: vmtkSeedSelector

Execute()[source]#
InitializeSeeds()[source]#
PickCallback(obj)[source]#
UndoCallback(obj)[source]#
class vampy.automatedPreprocessing.vmtk_pointselector.vmtkSeedSelector[source]#

Bases: object

Execute()[source]#
GetSurface()[source]#
GetTargetSeedIds()[source]#
SetSurface(surface)[source]#

Post-processing scripts#

vampy.automatedPostprocessing.compute_flow_and_simulation_metrics.compute_flow_and_simulation_metrics(folder, nu, dt, velocity_degree, T, times_to_average, save_frequency, start_cycle, step, average_over_cycles)[source]#

Computes several flow field characteristics and metrics for the velocity field. Reads in the .h5 soution stored in the ‘Solution’ folder within the results’ folder, and stores the post-processing files in the FlowMetrics folder.

Parameters
  • folder (str) – The folder path where the data files are located.

  • nu (float) – The kinematic viscosity value.

  • dt (float) – The time step size.

  • velocity_degree (int) – The degree of the velocity function space.

  • T (float) – Length of one cardiac cycle, in [ms]

  • times_to_average (list) – A list of time points to perform phase averaging. Needs to be in the inverval [0,T).

  • save_frequency (int) – The frequency of saving data during the simulation.

  • start_cycle (int) – The starting cycle number for averaging.

  • step (int) – The step size for extracting data.

  • average_over_cycles (bool) – A flag indicating whether to perform cycle averaging.

vampy.automatedPostprocessing.compute_flow_and_simulation_metrics.compute_u_avg(dataset_names, file_path_u_avg, file_u, n_cycles, saved_time_steps_per_cycle, start_cycle, u, u_avg)[source]#

Iterate over saved time steps and compute average velocity based on save sampling parameters

Parameters
  • dataset_names (dict) – Contains timestep and respective names of solution timesteps

  • file_path_u_avg (str) – File path to average velocity solution

  • file_u (str) – File path to velocity solution

  • n_cycles (int) – Number of cardiac cycles

  • saved_time_steps_per_cycle (int) – Determines number of saved time steps per cycle

  • start_cycle (int) – Determines which cardiac cycle to start sampling from

  • u (Function) – Function storing velocity

  • u_avg (Function) – Function for storing average velocity

  • start_cycle – Determines which cardiac cycle to start from for post-processing

vampy.automatedPostprocessing.compute_flow_and_simulation_metrics.define_functions_and_iterate_dataset(time_to_average, dataset, dataset_avg, dt, file_u, file_path_u_avg, file_path_u, mesh, nu, number_of_files, DG, V, Vv, cycles_to_average, saved_time_steps_per_cycle)[source]#

Defines functions and vector functions for all metrics to be computed, iterates through dataset and computes several flow and simulation metrics. After iteration, saves the metrics to .xdmf file format.

Parameters
  • time_to_average (int) – Time in ms to average over

  • number_of_files (int) – Number of dataset files to iterate over

  • DG (FunctionSpace) – Discontinous Galerkin function space

  • V (VectorFunctionSpace) – Continous Galerkin function space for vector solutions

  • Vv (FunctionSpace) – Continous Galerkin function space for scalar solutions

  • mesh (Mesh) – Volumetric mesh of computational domain

  • dataset (dict) – Contains velocity solution dict keys

  • dataset_avg (dict) – Contains average velocity solution dict keys

  • file_path_u (str) – File path to velocity solution

  • file_path_u_avg (str) – File path to average velocity solution

  • file_u (str) – File path to velocity solution

  • saved_time_steps_per_cycle (int) – Determines number of saved time steps per cycle

  • folder (str) – Path to simulation results

  • nu (float) – Viscosity

  • dt (float) – Time step in [ms]

  • cycles_to_average (list) – List of which cycles to average

vampy.automatedPostprocessing.compute_flow_and_simulation_metrics.get_files_for_cycle_averaging(dataset_names, file_path_u_avg, number_of_cycles, saved_time_steps_per_cycle, start_cycle)[source]#

Retrieves the dataset dictionaries for cycle averaging.

Parameters
  • dataset_names (list) – A list of dataset names.

  • file_path_u_avg (str) – The file path for averaging the dataset.

  • number_of_cycles (int) – The total number of cycles.

  • saved_time_steps_per_cycle (int) – The number of time steps saved per cycle.

  • start_cycle (int) – The starting cycle number for averaging.

Returns

A dictionary containing dataset names for cycle averaging. dataset_dict_avg (dict): A dictionary containing averaged dataset names. number_of_files (int): The total number of files.

Return type

dataset_dict (dict)

vampy.automatedPostprocessing.compute_flow_and_simulation_metrics.get_files_for_phase_averaging(times_to_average, dt, save_frequency, step, number_of_cycles, start_cycle, saved_time_steps_per_cycle, dataset_names)[source]#

Generates dataset dictionaries for phase averaging based on the selected times.

Parameters
  • times_to_average (list) – A list of time points to perform phase averaging.

  • dt (float) – The time step size.

  • save_frequency (int) – The frequency of saving data during the simulation.

  • step (int) – The step size for extracting data.

  • number_of_cycles (int) – The total number of cycles.

  • start_cycle (int) – The starting cycle number for averaging.

  • saved_time_steps_per_cycle (int) – The number of time steps saved per cycle.

  • dataset_names (list) – A list of dataset names.

Returns

A dictionary mapping time points to corresponding dataset names. dataset_dict_avg (dict): A dictionary mapping time points to averaged dataset names. number_of_files (int): The total number of files.

Return type

dataset_dict (dict)

vampy.automatedPostprocessing.compute_flow_and_simulation_metrics.main_metrics()[source]#
vampy.automatedPostprocessing.compute_hemodynamic_indices.compute_hemodynamic_indices(folder, nu, rho, dt, T, velocity_degree, save_frequency, start_cycle, step, average_over_cycles)[source]#

Loads velocity fields from completed CFD simulation, and computes and saves the following hemodynamic quantities: (1) WSS - Wall shear stress (2) TAWSS - Time averaged wall shear stress (3) TWSSG - Temporal wall shear stress gradient (4) OSI - Oscillatory shear index (5) RRT - Relative residence time (6) ECAP - Endothelial cell activation potential

The resulting wall shear stress will be in units Pascal [Pa], given that the provided density (rho) is in [kg/m^3], the time step (dt) is in [ms], and viscosity (nu) is in [mm^2/ms].

Parameters
  • velocity_degree (int) – Finite element degree of velocity

  • folder (str) – Path to results from simulation

  • nu (float) – Kinematic viscosity

  • rho (float) – Fluid density

  • dt (float) – Time step of simulation

  • T (float) – Length of one cardiac cycle, in [ms]

  • save_frequency (int) – Frequency that velocity has been stored

  • start_cycle (int) – Determines which cardiac cycle to start from for post-processing

  • step (int) – Step size determining number of times data is sampled

  • average_over_cycles (bool) – Averages over cardiac cycles if True

vampy.automatedPostprocessing.compute_hemodynamic_indices.main_hemo()[source]#
vampy.automatedPostprocessing.compute_velocity_and_pressure.compute_velocity_and_pressure(folder, dt, save_frequency, velocity_degree, pressure_degree, step)[source]#

Loads velocity and pressure from compressed .h5 CFD solution and converts and saves to .xdmf format for visualization (in e.g. ParaView).

Parameters
  • folder (str) – Path to results from simulation

  • dt (float) – Time step of simulation

  • save_frequency (int) – Frequency that velocity and pressure has been stored

  • velocity_degree (int) – Finite element degree of velocity

  • pressure_degree (int) – Finite element degree of pressure

  • step (int) – Step size determining number of times data is sampled

vampy.automatedPostprocessing.compute_velocity_and_pressure.main_convert()[source]#
class vampy.automatedPostprocessing.postprocessing_common.STRESS(u, p, nu, mesh)[source]#

Bases: object

Computes the stress on a given mesh based on provided velocity and pressure fields. Note that the pressure term is unused since it disappears in the process of extracting the tangential component.

FIXME: Currently works for P1P1, but not for higher order finite elements (e.g. P2P1)

norm_l2(u)[source]#

Compute norm of vector u in expression form :param u: Function to compute norm of :type u: Function

Returns

Norm as expression

Return type

norm (Power)

vampy.automatedPostprocessing.postprocessing_common.epsilon(u)[source]#

Computes the strain-rate tensor :param u: Velocity field :type u: Function

Returns

Strain rate tensor of u

Return type

epsilon (Function)

vampy.automatedPostprocessing.postprocessing_common.get_dataset_names(data_file, num_files=100000, step=1, start=0, print_info=True, vector_filename='/velocity/vector_%d')[source]#

Read velocity fields datasets and extract names of files

Parameters
  • data_file (HDF5File) – File object of velocity

  • num_files (int) – Number of files

  • step (int) – Step between each data dump

  • start (int) – Step to start on

  • print_info (bool) – Prints info about data if true

  • vector_filename (str) – Name of velocity files

Returns

List of data file names

Return type

names (list)

vampy.automatedPostprocessing.postprocessing_common.rate_of_dissipation(dissipation, u, v, mesh, h, nu)[source]#

Computes rate of dissipation

Parameters
  • dissipation (Function) – Function to save rate of dissipation to

  • u (Function) – Function for velocity field

  • v (Function) – Test function for velocity

  • mesh – Mesh to compute dissipation rate on

  • h (float) – Cell diameter of mesh

  • nu (float) – Viscosity

vampy.automatedPostprocessing.postprocessing_common.rate_of_strain(strain, u, v, mesh, h)[source]#

Computes rate of strain

Parameters
  • strain (Function) – Function to save rate of strain to

  • u (Function) – Function for velocity field

  • v (Function) – Test function for velocity

  • mesh – Mesh to compute strain rate on

  • h (float) – Cell diameter of mesh

vampy.automatedPostprocessing.postprocessing_common.read_command_line()[source]#

Read arguments from commandline

vampy.automatedPostprocessing.visualize_probes.main_probe()[source]#
vampy.automatedPostprocessing.visualize_probes.visualize_probes(case_path, probe_saving_frequency, show_figure=True, save_figure=False)[source]#

Loads probe points from completed CFD simulation, and visualizes velocity (magnitude) and pressure at respective probes Assuming results are in mm/s.

Parameters
  • probe_saving_frequency (int) – Interval between saving probes.

  • case_path (str) – Path to results from simulation

  • show_figure (bool) – Shows figure if True

  • save_figure (bool) – Saves figure if True

CFD simulation scripts#

vampy.simulation.Artery.beta(err, p)[source]#

Adjusted choice of beta for the dual-pressure boundary condition. Ramped up to desired value if flow rate error (err) increases

Parameters
  • err (float) – Flow split error

  • p (float) – Pressure value

Returns

Variable factor in flow split method

Return type

beta (float)

vampy.simulation.Artery.create_bcs(t, NS_expressions, V, Q, area_ratio, area_inlet, mesh, mesh_path, nu, id_in, id_out, pressure_degree, **NS_namespace)[source]#
vampy.simulation.Artery.mesh(mesh_path, **NS_namespace)[source]#
vampy.simulation.Artery.pre_solve_hook(mesh, V, Q, newfolder, mesh_path, restart_folder, velocity_degree, cardiac_cycle, save_solution_after_cycle, dt, **NS_namespace)[source]#
vampy.simulation.Artery.problem_parameters(commandline_kwargs, NS_parameters, NS_expressions, **NS_namespace)[source]#

Problem file for running CFD simulation in arterial models consisting of one inlet, and two or more outlets. A Womersley velocity profile is applied at the inlet, and a flow split pressure condition is applied at the outlets, following [1]. Flow rate for the inlet condition, and flow split values for the outlets are computed from the pre-processing script automatedPreProcessing.py. The simulation is run for two cycles (adjustable), but only the results/solutions from the second cycle are stored to avoid non-physiological effects from the first cycle. One cardiac cycle is set to 0.951 s from [2], and scaled by a factor of 1000, hence all parameters are in [mm] or [ms].

[1] Gin, Ron, Anthony G. Straatman, and David A. Steinman. “A dual-pressure boundary condition for use in

simulations of bifurcating conduits.” J. Biomech. Eng. 124.5 (2002): 617-619.

[2] Hoi, Yiemeng, et al. “Characterization of volumetric flow rate waveforms at the carotid bifurcations of older

adults.” Physiological measurement 31.3 (2010): 291.

vampy.simulation.Artery.temporal_hook(u_, p_, mesh, tstep, dump_probe_frequency, eval_dict, newfolder, id_in, id_out, boundary, n, save_solution_frequency, NS_parameters, NS_expressions, area_ratio, t, save_solution_at_tstep, U, area_inlet, nu, u_mean0, u_mean1, u_mean2, **NS_namespace)[source]#
vampy.simulation.Artery.theend_hook(u_mean, u_mean0, u_mean1, u_mean2, T, dt, save_solution_at_tstep, save_solution_frequency, **NS_namespace)[source]#
vampy.simulation.Artery.update_pressure_condition(NS_expressions, area_ratio, boundary, id_in, id_out, mesh, n, tstep, u_)[source]#

Use a dual-pressure boundary condition as pressure condition at outlet.

vampy.simulation.Atrium.create_bcs(NS_expressions, mesh, mesh_path, nu, t, V, Q, id_in, id_out, **NS_namespace)[source]#
vampy.simulation.Atrium.mesh(mesh_path, **NS_namespace)[source]#
vampy.simulation.Atrium.pre_solve_hook(V, Q, cardiac_cycle, dt, save_solution_after_cycle, mesh_path, mesh, newfolder, velocity_degree, restart_folder, **NS_namespace)[source]#
vampy.simulation.Atrium.problem_parameters(commandline_kwargs, NS_parameters, scalar_components, Schmidt, NS_expressions, **NS_namespace)[source]#

Problem file for running CFD simulation in left atrial models consisting of arbitrary number of pulmonary veins (PV) (normally 3 to 7), and one outlet being the mitral valve. A Womersley velocity profile is applied at the inlets, where the total flow rate is split between the area ratio of the PVs. The mitral valve is considered open with a constant pressure of p=0. Flow rate and flow split values for the inlet condition are computed from the pre-processing script automatedPreProcessing.py. The simulation is run for two cycles (adjustable), but only the results/solutions from the second cycle are stored to avoid non-physiological effects from the first cycle. One cardiac cycle is set to 0.951 s from [1], and scaled by a factor of 1000, hence all parameters are in [mm] or [ms].

[1] Hoi, Yiemeng, et al. “Characterization of volumetric flow rate waveforms at the carotid bifurcations of older

adults.” Physiological measurement 31.3 (2010): 291.

vampy.simulation.Atrium.temporal_hook(mesh, dt, t, save_solution_frequency, u_, NS_expressions, id_in, tstep, newfolder, eval_dict, dump_probe_frequency, p_, save_solution_at_tstep, nu, D_mitral, U, u_mean0, u_mean1, u_mean2, **NS_namespace)[source]#
vampy.simulation.Atrium.theend_hook(u_mean, u_mean0, u_mean1, u_mean2, T, dt, save_solution_at_tstep, save_solution_frequency, **NS_namespace)[source]#
class vampy.simulation.Womersley.WomersleyComponent(radius, center, normal, normal_component, period, nu, element, Q=None, V=None)[source]#

Bases: UserExpression

eval(value, x)[source]#
set_t(t)[source]#
vampy.simulation.Womersley.compute_boundary_geometry_acrn(mesh, ind, facet_domains)[source]#
vampy.simulation.Womersley.fourier_coefficients(x, y, T, N)[source]#

From x-array and y-spline and period T, calculate N complex Fourier coefficients.

vampy.simulation.Womersley.make_womersley_bcs(t, Q, nu, center, radius, normal, element, coeffstype='Q', N=1001, num_fourier_coefficients=20, Cn=None, **NS_namespace)[source]#

Generate a list of expressions for the components of a Womersley profile. Users can specify either the flow rate or fourier coefficients of the flow rate depending on if Cn is None or not.

vampy.simulation.Womersley.x_to_r2(x, c, n)[source]#

Compute r**2 from a coordinate x, center point c, and normal vector n.

r is defined as the distance from c to x’, where x’ is the projection of x onto the plane defined by c and n.

vampy.simulation.simulation_common.dump_probes(eval_dict, newfolder, tstep)[source]#

Save variables along the centerline for CFD simulation diagnostics and light-weight post-processing.

Parameters
  • eval_dict (dict) – Dictionary with probe arrays for velocity components and pressure along the centerline.

  • newfolder (str) – The path to the folder where the probe data will be saved.

  • tstep (float) – The current time step.

Returns

None

vampy.simulation.simulation_common.get_file_paths(folder, additional_variables=[])[source]#

Create folder where data and solutions (velocity, mesh, pressure) is stored

Parameters
  • folder (str) – Path to data storage location

  • additional_variables (list) – List of additional variables to store

Returns

Contains filepaths for respective solution files

Return type

files (dict)

vampy.simulation.simulation_common.print_mesh_information(mesh)[source]#

Print geometric information about the volumetric mesh

Parameters

mesh (dolfin.cpp.mesh.Mesh) – Volumetric mesh

vampy.simulation.simulation_common.store_u_mean(T, dt, save_solution_at_tstep, save_solution_frequency, u_mean, u_mean0, u_mean1, u_mean2, NS_parameters)[source]#

Store the time averaged velocity into file u_mean

Parameters
  • T (float) – End time

  • dt (float) – Time step size

  • save_solution_at_tstep (int) – Time step when solution is to be saved

  • save_solution_frequency (int) – Frequency of solution saving

  • u_mean (Function) – Vector function for storing vector solution of average velocity

  • u_mean0 (Function) – Function storing x-component

  • u_mean1 (Function) – Function storing y-component

  • u_mean2 (Function) – Function storing z-component

vampy.simulation.simulation_common.store_velocity_and_pressure_h5(NS_parameters, U, p_, tstep, u_, u_mean0, u_mean1, u_mean2, D=None, du_=None)[source]#

Store the velocity and pressure values to an HDF5 file.

Parameters
  • NS_parameters (dict) – A dictionary containing the parameters for Navier-Stokes equations.

  • U (Function) – A vector function space to assign the velocity components.

  • p (Function) – The pressure function.

  • tstep (int) – The current time step.

  • u (List[Function]) – A list containing the velocity components.

  • u_mean0 (Function) – The accumulated x-component of the velocity.

  • u_mean1 (Function) – The accumulated y-component of the velocity.

  • u_mean2 (Function) – The accumulated z-component of the velocity.

  • du (List[Function]) – List of deformation components

  • D (Function) – A vector function space to assign the deformation components.

Returns

None