purkinje_uv package¶
Submodules¶
purkinje_uv.branch module¶
Module defining the Branch class for fractal tree growth on a mesh.
This module contains the Branch class, which represents a single branch in the fractal tree.
- class purkinje_uv.branch.Branch(mesh, init_node, init_dir, init_tri, length, angle, w, nodes, brother_nodes, Nsegments)
Bases:
object
Represents a branch of the fractal tree on a mesh.
- Parameters:
mesh (Mesh)
init_node (int)
init_dir (NDArray[Any])
init_tri (int)
length (float)
angle (float)
w (float)
nodes (Nodes)
brother_nodes (Sequence[int])
Nsegments (int)
- mesh
Mesh on which the branch grows.
- Type:
Mesh
- queue
Coordinates of points queued for growth.
- Type:
list[NDArray[Any]]
- triangles
Triangle indices for each queued point.
- Type:
list[int]
- nodes
Node indices in the global Nodes manager.
- Type:
list[int]
- growing
Whether the branch is still growing.
- Type:
bool
- add_node_to_queue(mesh, init_node, dir_vec)
Project a new node onto the mesh and add it to the growth queue.
This method projects a direction from a starting point onto the mesh surface. If the projected point lies within a mesh triangle, it is appended to the queue and triangles list.
- Parameters:
mesh (Mesh) – Mesh on which the branch grows.
init_node (NDArray[Any]) – Coordinates of the last node in the branch.
dir_vec (NDArray[Any]) – Direction vector from the initial node to the new node.
- Returns:
True if the projected node lies within a mesh triangle; False otherwise.
- Return type:
bool
purkinje_uv.edge module¶
Module defining the Edge class for graph edges in a Purkinje tree.
This module provides the Edge class, representing a connection between two nodes, including its geometric direction and optional parent/branch relationships.
- class purkinje_uv.edge.Edge(n1, n2, nodes, parent, branch)
Bases:
object
Represents an edge between two nodes in a graph.
Each Edge connects node n1 to node n2, computes the normalized direction vector between them, and optionally tracks a parent edge and branch association.
- Parameters:
n1 (int)
n2 (int)
nodes (Sequence[NDArray[Any]])
parent (Optional[int])
branch (Optional[int])
- n1
The ID of the first node.
- Type:
int
- n2
The ID of the second node.
- Type:
int
- dir
Normalized direction vector from n1 to n2.
- Type:
NDArray[Any]
- parent
The parent edge ID, if any.
- Type:
Optional[int]
- branch
The branch ID, if any.
- Type:
Optional[int]
purkinje_uv.fractal_tree_parameters module¶
Module defining the FractalTreeParameters class for configuring fractal tree generation.
This module provides the FractalTreeParameters dataclass, which holds all settings for generating a fractal tree structure.
- class purkinje_uv.fractal_tree_parameters.FractalTreeParameters(meshfile=None, init_node_id=0, second_node_id=1, init_length=0.1, N_it=10, length=0.1, branch_angle=0.15, w=0.1, l_segment=0.01, fascicles_angles=<factory>, fascicles_length=<factory>)
Bases:
object
Holds settings for generating a fractal tree structure.
- Parameters:
meshfile (str | None)
init_node_id (int)
second_node_id (int)
init_length (float)
N_it (int)
length (float)
branch_angle (float)
w (float)
l_segment (float)
fascicles_angles (List[float])
fascicles_length (List[float])
- meshfile
Path to the mesh file (e.g., OBJ); if
None
, no mesh is preloaded.- Type:
str | None
- init_node_id
Index of the initial node in the mesh.
- Type:
int
- second_node_id
Index of the second node, which sets the initial growth direction.
- Type:
int
- init_length
Length of the first branch.
- Type:
float
- N_it
Number of branch generations (iterations).
- Type:
int
- length
Mean length of tree branches.
- Type:
float
- branch_angle
Angle (in radians) between consecutive branches.
- Type:
float
- w
Weight parameter for branch divergence.
- Type:
float
- l_segment
Approximate length of each branch segment.
- Type:
float
- fascicles_angles
Angles (in radians) for each fascicle branch.
- Type:
List[float]
- fascicles_length
Lengths for each fascicle branch; matches
fascicles_angles
.- Type:
List[float]
- N_it: int
- branch_angle: float
- fascicles_angles: List[float]
- fascicles_length: List[float]
- classmethod from_dict(data)
Construct parameters from a dictionary.
Unknown keys are ignored (emits a warning). Validation occurs in
__post_init__
.- Parameters:
data (Mapping[str, Any])
- Return type:
FractalTreeParameters
- classmethod from_json(s)
Construct parameters from a JSON string.
- Parameters:
s (str)
- Return type:
FractalTreeParameters
- init_length: float
- init_node_id: int
- l_segment: float
- length: float
- meshfile: str | None
- second_node_id: int
- to_dict()
Return a JSON-serializable dictionary of the parameters.
- Return type:
Dict[str, Any]
- to_json(indent=2)
Serialize parameters to a JSON string.
- Parameters:
indent (int | None) – Indentation level for pretty printing. If
None
, the result is compact.- Return type:
str
- to_json_file(path, indent=2)
Write parameters to a JSON file at
path
.- Parameters:
path (str)
indent (int | None)
- Return type:
None
- w: float
purkinje_uv.fractal_tree module¶
Module defining the FractalTree class to generate fractal trees within a mesh domain.
This module implements UV-based fractal tree growth using geometric rules, collision detection, and iterative branching to create a tree structure embedded in a 3D mesh.
- class purkinje_uv.fractal_tree.FractalTree(params)
Bases:
object
Fractal-tree generator on a UV-mapped surface.
This class reproduces the legacy algorithm exactly (VTK-based point-in-mesh checks, UV-scaling step sizes, and queue behavior), keeping helpers local so the module is self-contained.
- Parameters:
params (FractalTreeParameters)
- grow_tree()
Generate the Purkinje fractal tree.
- Return type:
None
- save(filename)
Write the generated line mesh to a VTU file.
- Parameters:
filename (str)
- Return type:
None
purkinje_uv.mesh module¶
Module defining the Mesh class for 3D triangular surface meshes.
- This module provides:
Loading from OBJ or direct arrays.
Computation of normals, centroids, connectivity.
KD-tree spatial queries.
FEM routines (B-matrix, stiffness, mass, force).
Designed for fractal tree growth and surface-based FEM analysis.
- class purkinje_uv.mesh.Mesh(filename=None, verts=None, connectivity=None)
Bases:
object
Handle 3D triangular surface meshes.
Supports geometry, topology, and finite-element operations. It computes normals, centroids, boundary edges; builds KD-trees; and provides FEM element routines, geodesic/Laplace solvers, UV mapping, and interpolation utilities.
- Parameters:
filename (Optional[str]) – Path to OBJ file to load mesh from.
verts (Optional[NDArray[Any]]) – Vertex coordinates (n_nodes×3).
connectivity (Optional[NDArray[Any]]) – Triangle indices (n_triangles×3).
- verts
Vertex array, shape (n_nodes, 3).
- Type:
NDArray[Any]
- connectivity
Triangle indices, shape (n_triangles, 3).
- Type:
NDArray[Any]
- normals
Triangle normals, shape (n_triangles, 3).
- Type:
NDArray[Any]
- node_to_tri
Node→[triangle indices].
- Type:
DefaultDict[int, List[int]]
- tree
KD-tree over verts for nearest-node queries.
- Type:
cKDTree
- centroids
Triangle centroids, shape (n_triangles, 3).
- Type:
NDArray[Any]
- boundary_edges
Edges on mesh boundary.
- Type:
Optional[List[Tuple[int,int]]]
- uv
UV coordinates per node, shape (n_nodes, 2).
- Type:
Optional[NDArray[Any]]
- triareas
Triangle areas, shape (n_triangles,).
- Type:
Optional[NDArray[Any]]
- uvscaling
UV scaling metric per triangle.
- Type:
Optional[NDArray[Any]]
- Bmatrix(element)
Compute the B-matrix and Jacobian determinant for a triangle.
- Parameters:
element (int) – Triangle index.
- Returns:
B (2×3 array): Strain-displacement matrix.
J (float): Twice the triangle area (Jacobian determinant).
- Return type:
Tuple[NDArray[Any], float]
- ForceVector(B, J, X)
Compute the local force vector for a triangle.
- Parameters:
B (NDArray[Any]) – B-matrix from Bmatrix.
J (float) – Jacobian determinant.
X (NDArray[Any]) – Gradient vector from gradient.
- Returns:
Length-3 force vector.
- Return type:
NDArray[Any]
- MassMatrix(J)
Compute the 3×3 (consistent) local mass matrix for a triangle.
- Parameters:
J (float) – Jacobian determinant.
- Returns:
3×3 mass matrix.
- Return type:
NDArray[Any]
- StiffnessMatrix(B, J)
Compute the local stiffness matrix for a triangle.
- Parameters:
B (NDArray[Any]) – B-matrix from Bmatrix.
J (float) – Jacobian determinant.
- Returns:
3×3 stiffness matrix.
- Return type:
NDArray[Any]
- boundary_edges: List[Tuple[int, int]] | None
- centroids: ndarray[tuple[int, ...], dtype[Any]]
- computeGeodesic(nodes, nodeVals, filename=None, K=None, M=None, dt=10.0)
Compute geodesic distances using the heat method (FEM).
- Parameters:
nodes (Sequence[int]) – Indices of fixed-temperature nodes.
nodeVals (Sequence[float]) – Temperature values at nodes.
filename (Optional[str]) – VTU output path.
K (Optional[sp.spmatrix]) – Preassembled stiffness matrix.
M (Optional[sp.spmatrix]) – Preassembled mass matrix.
dt (float) – Time-step for the heat diffusion.
- Returns:
ATglobal: Per-node geodesic distance.
Xs: Per-triangle gradient directions.
- Return type:
Tuple[NDArray[Any], NDArray[Any]]
- computeLaplace(nodes, nodeVals, filename=None)
Solve Laplace’s equation with Dirichlet boundary conditions.
- Parameters:
nodes (Sequence[int]) – Dirichlet node indices.
nodeVals (Sequence[float] | NDArray[Any]) – Boundary values.
filename (Optional[str]) – VTU output path.
- Returns:
Solution vector of length n_nodes.
- Return type:
NDArray[Any]
- computeLaplacian()
Assemble global stiffness (K) and mass (M) matrices as CSR.
- Returns:
(K, M) in CSR format.
- Return type:
Tuple[sp.spmatrix, sp.spmatrix]
- compute_triareas()
Compute triangle areas (J/2) and store in self.triareas.
- Return type:
None
- compute_uvscaling()
Compute and validate UV-scaling metric per triangle.
- Return type:
None
- connectivity: ndarray[tuple[int, ...], dtype[Any]]
- detect_boundary()
Identify boundary edges (edges in exactly one triangle).
- Return type:
None
- gradient(element, u)
Compute the gradient of a scalar field over a triangle.
- Parameters:
element (int) – Triangle index.
u (NDArray[Any]) – Field values at the triangle’s three vertices.
- Returns:
3-vector of gradients in 3D space.
- Return type:
NDArray[Any]
- loadOBJ(filename)
Read a Wavefront .obj mesh file and return (verts, connectivity).
- Parameters:
filename (str) – Path to the .obj file.
- Returns:
verts: Array of shape (n_vertices, 3)
connectivity: Array of shape (n_triangles, 3)
- Return type:
Tuple[NDArray[Any], NDArray[Any]]
- node_to_tri: DefaultDict[int, List[int]]
- normals: ndarray[tuple[int, ...], dtype[Any]]
- project_new_point(point, verts_to_search=1)
Project a point onto the mesh and find its containing triangle.
- Parameters:
point (NDArray[Any]) – Coordinates to project.
verts_to_search (int) – Number of nearby vertices to search.
- Returns:
Projected point, triangle index (−1 if outside), and barycentric coords (r, t).
- Return type:
Tuple[NDArray[Any], int, float, float]
- project_point_check(point, node)
This function projects any point to the surface defined by the mesh.
- Parameters:
point (array) – coordinates of the point to project.
node (int) – index of the most close node to the point
- Returns:
the coordinates of the projected point that lies in the surface. intriangle (int): the index of the triangle where the projected point lies. If the point is outside surface, intriangle=-1.
- Return type:
projected_point (array)
- tree: cKDTree
- tri2node_interpolation(cell_field)
Interpolate triangle-based field to nodes by area-weighting.
- Parameters:
cell_field (NDArray[Any]) – Per-triangle values.
- Returns:
Per-node interpolated values.
- Return type:
List[float]
- triareas: ndarray[tuple[int, ...], dtype[Any]] | None
- uv: ndarray[tuple[int, ...], dtype[Any]] | None
- uv_bc()
Generate UV boundary loop and boundary conditions.
- Returns:
around_nodes, bc_u, bc_v.
- Return type:
Tuple[List[int], NDArray[Any], NDArray[Any]]
- uvmap(filename=None)
Compute UV coordinates by solving Laplace’s equation.
- Parameters:
filename (Optional[str]) – VTU export path for u and v fields.
- Return type:
None
- uvscaling: ndarray[tuple[int, ...], dtype[Any]] | None
- verts: ndarray[tuple[int, ...], dtype[Any]]
- writeVTU(filename, point_data=None, cell_data=None)
Export this mesh (and optional point/cell data) in VTU format.
Accepts NumPy or CuPy arrays; all data are converted to NumPy on CPU before writing. cell_data is normalized to meshio’s expected shape (dict of name -> list-of-arrays, per cell block).
- Parameters:
filename (str) – Output path (e.g.,
"mesh.vtu"
).point_data (Dict[str, ndarray[tuple[int, ...], dtype[Any]]] | None) – Optional dict of per-node arrays (shape (n_nodes,) or (n_nodes, k)).
cell_data (Dict[str, ndarray[tuple[int, ...], dtype[Any]]] | None) – Optional dict of per-triangle arrays (shape (n_tris,) or (n_tris, k)).
- Raises:
ValueError – If provided data have incompatible lengths.
Exception – If the underlying mesh writer fails.
- Return type:
None
purkinje_uv.nodes module¶
Module defining the Nodes class for managing tree nodes and spatial queries.
This module provides the Nodes class, which stores branch nodes and offers methods to compute distances, collisions, and gradients via KD‐trees.
- class purkinje_uv.nodes.Nodes(init_node)
Bases:
object
Manage nodes and compute spatial queries for a tree structure.
The Nodes class stores branch nodes and provides methods to compute distances, collisions, and gradients using k-d trees.
- Parameters:
init_node (NDArray[Any])
- nodes
Coordinates of all nodes.
- Type:
List[NDArray[Any]]
- last_node
Index of the most recently added node.
- Type:
int
- end_nodes
Indices of terminal nodes (not connected).
- Type:
List[int]
- tree
KD-tree of all nodes for nearest-neighbor queries.
- Type:
cKDTree
- collision_tree
KD-tree excluding certain nodes for collision checks.
- Type:
cKDTree
- add_nodes(queue)
Append a sequence of new nodes and rebuild the KD-tree.
- Parameters:
queue (Sequence[NDArray[Any]]) – Coordinates of nodes to add.
- Returns:
Indices of the newly added nodes.
- Return type:
List[int]
- collision(point)
Find the nearest node (excluding excluded ones) to a query point.
- Parameters:
point (Union[NDArray[Any], Sequence[float]]) – Query coordinates.
- Returns:
(node_index, distance) to the closest node.
- Return type:
Tuple[int, float]
- collision_tree: cKDTree
- distance_from_node(node)
Compute distance from one node to its nearest neighbor in the tree.
- Parameters:
node (int) – Index of the node to query.
- Returns:
Distance to the closest other node.
- Return type:
float
- distance_from_point(point)
Compute distance from an arbitrary point to the nearest node.
- Parameters:
point (Union[NDArray[Any], Sequence[float]]) – Query coordinates.
- Returns:
Distance to the closest node.
- Return type:
float
- end_nodes: List[int]
- gradient(point)
Approximate the gradient of the distance field at a point.
Uses a central difference if needed, but by default returns a unit vector pointing away from the nearest node.
- Parameters:
point (Union[NDArray[Any], Sequence[float]]) – Query coordinates.
- Returns:
(dx, dy, dz) gradient components.
- Return type:
NDArray[Any]
- nodes: List[ndarray[tuple[int, ...], dtype[Any]]]
- nodes_to_consider_keys: list[int]
- tree: cKDTree
- update_collision_tree(nodes_to_exclude)
Rebuild the collision KD-tree excluding specified nodes.
If all nodes are excluded, inserts a distant dummy point so the KD-tree remains non-empty. Returns
True
on success; on failure, preserves the previous collision tree and returnsFalse
.- Parameters:
nodes_to_exclude (Sequence[int]) – Indices to omit from collision checks.
- Returns:
True
if the collision KD-tree was rebuilt and installed,False
if an error occurred (previous tree remains active).- Return type:
bool
purkinje_uv.purkinje_tree module¶
Module for eikonal activation on a Purkinje fiber network.
This module defines the PurkinjeTree class, which wraps a line-based Purkinje fiber mesh, computes activation times using the Fast Iterative Method (FIM), and provides utilities for exporting results to VTK and meshio formats.
- class purkinje_uv.purkinje_tree.PurkinjeTree(nodes, connectivity, end_nodes)
Bases:
object
Perform eikonal activation on a Purkinje network using FIM solver.
- Parameters:
nodes (Sequence[ndarray[tuple[int, ...], dtype[Any]]])
connectivity (Sequence[Sequence[int]])
end_nodes (Sequence[int])
- activate_fim(x0, x0_vals, return_only_pmj=True)
Compute activation times using the Fast Iterative Method (FIM).
- Parameters:
x0 (NDArray[Any]) – Array of source node indices.
x0_vals (NDArray[Any]) – Activation values at each source node.
return_only_pmj (bool) – If True, return activation times only at PMJ nodes.
- Returns:
Activation times for all nodes or only PMJ nodes.
- Return type:
NDArray[Any]
- extract_edges()
Return the edge connectivity array for the Purkinje tree.
- Returns:
Array of shape (n_edges, 2) listing node index pairs.
- Return type:
NDArray[Any]
- extract_pmj_counter()
Compute PMJ node indices using a pure Python counter on connectivity.
- Returns:
Nodes with degree equal to one.
- Return type:
List[int]
- extract_pmj_np_bincount()
Compute PMJ node indices by numpy bincount on flattened connectivity.
- Returns:
Nodes with degree equal to one.
- Return type:
NDArray[Any]
- extract_pmj_np_unique()
Compute PMJ node indices via numpy unique and count filtering.
- Returns:
Nodes with degree equal to one.
- Return type:
NDArray[Any]
- get_pmjs_activation()
Retrieve activation times at the Purkinje-myocardial junction nodes.
- Returns:
Activation times indexed by PMJ node order.
- Return type:
NDArray[Any]
- save(fname)
Write the current activation state to a VTK UnstructuredGrid file.
- Parameters:
fname (str) – Path to the output .vtu file.
- Return type:
None
- save_meshio(fname, point_data=None, cell_data=None)
Export the Purkinje tree as a meshio Mesh with line cells.
- Parameters:
fname (str)
point_data (Dict[str, Any] | None)
cell_data (Dict[str, Any] | None)
- Return type:
None
- save_pmjs(fname)
Export PMJ nodes and their activation values as a VTP file.
- Parameters:
fname (str) – Path to the output .vtp file.
- Return type:
None