Developer documentation
The developer documentation is organized by modules. This reflects the actual structure of the code. It further lists all private members of classes. This is useful for developers, but not for users.
In contrast, the user documentation directly lists all members of the package that the user can import directly from the hemtracer package. These members are defined in __init__.py
For example, the user imports the hemtracer.eulerian_flow_field.EulerianFlowField
class from the hemtracer
package as from hemtracer import EulerianFlowField
, not from the hemtracer.eulerian_flow_field
module. This should make working with the package easier for the user.
eulerian_flow_field
- class hemtracer.eulerian_flow_field.EulerianFlowField(vtk_filename)
Bases:
object
A class that represents a steady flow field in Eulerian coordinates with various field data. Unsteady flow fields are currently only supported in the form of MRF data. The mesh has to be unstructured and non-moving and the flow field has to be defined as node-centered velocity field (point data in VTK). The flow field is represented by a VTK file. It has to contain a velocity field. Additional field data, e.g., an Eulerian hemolysis solution, can be imported from other VTK files using the import_field_data method. The class provides functionality to interpolate field data to pathlines.
Initialize the EulerianField object with a VTK flow field name.
- Parameters:
vtk_filename (str) – Name of the VTK file containing the flow field.
- _compute_gradients_cell()
Compute the velocity gradient tensor as cell-centered field.
- Return type:
None
- _compute_gradients_point()
Compute the velocity gradient tensor as node-centered field.
- Return type:
None
- compute_distance_rotation(omega, x0=array([0., 0., 0.]))
Compute orthogonal distance from rotation axis. This is used for the transformation of the velocity field in a rotating frame of reference and for the tank-treading morphology model with rotation correction.
- Parameters:
omega (Vector3) – The angular velocity vector of the rotating frame of reference.
x0 (Vector3) – Some point on the rotation axis. Defaults to (0,0,0).
- Return type:
None
- compute_gradients(node_centered=False)
Compute the velocity gradient tensor (if it has not been computed before). The velocity gradient tensor can be computed as cell-centered field or node-centered field. The node-centered field is more expensive to compute, but it can be interpolated to pathlines automatically. The cell-centered field is cheaper to compute, but it has to be interpolated manually to pathlines.
- Parameters:
node_centered (bool) – If True, compute the gradient as node-centered field, else compute as cell-centered. Defaults to False.
- Return type:
None
- get_name_advection_velocity()
Get the name of the velocity field relevant for advection. In case of a rotating frame of reference, this is the MRF velocity field. Otherwise, this is the velocity field.
- Returns:
The velocity field name.
- Return type:
str
- get_name_distance_center()
Get the name of the orthogonal distance from the rotation axis. If no rotating frame of reference has been defined, this is None.
- Returns:
The orthogonal distance from the rotation axis.
- Return type:
str
- get_name_omega_frame()
Get the name of the angular velocity of the rotating frame of reference. If no rotating frame of reference has been defined, this is None.
- Returns:
The angular velocity of the rotating frame of reference.
- Return type:
str
- get_name_velocity_gradient()
Get the name of the velocity gradient field. If no velocity gradient has been computed, this is None.
- Returns:
The velocity gradient field name.
- Return type:
str
- get_vtk_flow_field()
Get the VTK flow field.
- Returns:
The VTK flow field.
- Return type:
vtkUnstructuredGrid
- import_field_data(vtk_file_name, field_names_point=None, field_names_cell=None)
Import additional field data from a VTK file (to be interpolated to the streamlines). The field data has to be defined on the same mesh as the flow field.
- Parameters:
vtk_file_name (str) – The name of the VTK file containing the additional field data.
field_names_point (List[str]) – A list of strings containing the names of the node-centered fields to import. If None, all fields are imported. If you do not want to import any node-centered fields, set field_names_point=[].
field_names_cell (List[str]) – A list of strings containing the names of the cell-centered fields to import. If None, all fields are imported. If you do not want to import any cell-centered fields, set field_names_cell=[].
- Return type:
None
- probe_field(x, field_names)
Extracts field information at specified points x, e.g., a pathline.
- Parameters:
x (List[Vector3]) – Positions in field at which to extract information.
field_names (List[str]) – List of names of field in data to extract.
- Returns:
A dictionary where keys are field names and values are the corresponding field information.
- Return type:
Dict[str, NDArray]
- reduce_to_relevant_fields(relevant_field_names)
Specify the names of the relevant fields in the VTK file. In particular, this could be velocity and all fields you want interpolated to the pathlines, e.g., some Eulerian solution for effective shear rate or local fluid shear. All other fields are removed. This can reduce memory requirements and speed up subsequent computations.
- Parameters:
relevant_field_names (List[str]) – A list of strings containing the names of the relevant fields.
- Return type:
None
- set_velocity_gradient_name(velocity_gradient_name)
Specify the name of the velocity gradient field in the VTK file, if available. The velocity gradient field can be defined as node-centered field (point data in VTK) or cell-centered field (cell data in VTK).
- Parameters:
velocity_gradient_name (str) – The name of the velocity gradient field in the VTK file.
- Return type:
None
- set_velocity_name(vtk_velocity_name)
Specify the name of the velocity field in the VTK file (velocity transformed to absolute frame for MRF fields). The velocity field has to be defined as node-centered field (point data in VTK)
- Parameters:
vtk_velocity_name (str) – The name of the velocity field in the VTK file.
- Return type:
None
- smooth_node_data(attribute_name)
Smooth certain node-centered data in the VTK file.
- Parameters:
attribute_name (str) – The name of the node-centered attribute to smooth.
- Return type:
None
- transform_mrf_data(rf_name, rf_rot, omega, x0=array([0., 0., 0.]))
Transform the velocity field of a rotating frame of reference.
- Parameters:
rf_name (str) – The name of the cell-centered (!) attribute in the vtk file that specifies the reference frame ID.
rf_rot (int) – The reference frame ID of the rotating frame of reference.
omega (Vector3) – The angular velocity vector of the rotating frame of reference.
x0 (Vector3) – Some point on the rotation axis. Defaults to (0,0,0).
- Return type:
None
pathlines
- class hemtracer.pathlines.Pathline(t, x, reason_for_termination='misc')
Bases:
object
Class representing a single pathline as a collection of PathlineAttribute objects. Every pathline has an attribute ‘Position’, which is stored by default as definition of the pathline. Additional attributes can be added using the add_attribute function.
Constructing an object creates its ‘Position’ attribute and stores it in attributes list.
- Parameters:
t (List[float]) – The integration times along the pathline.
x (List[Vector3]) – The positions along the pathline.
reason_for_termination (str) – The reason for termination of the pathline integration, defaults to ‘misc’.
- _attributes
A list of PathlineAttribute objects.
- _reason_for_termination
The reason for termination of the pathline integration. Can be ‘out of domain’, ‘not initialized’, ‘unexpected value’, ‘out of length’, ‘out of steps’, ‘stagnation’, ‘misc’.
- _t0
The initial integration time of the pathline.
- _tend
The final integration time of the pathline.
- add_attribute(t, val, name, interpolation_scheme='linear')
Adds an attribute to the pathline.
- Parameters:
t (ArrayLike) – The integration times along the pathline. If not in the range [t0, tend], the integration times are clipped.
val (ArrayLike) – The values of the attribute along the pathline (first axis must correspond to integration times). If the attribute is vector-valued, the second axis must correspond to the components of the vector.
name (str) – The name of the attribute.
interpolation_scheme (str) – The interpolation scheme to use along the pathline. Defaults to ‘linear’.
- Return type:
None
- get_attribute(name)
Returns the attribute with the specified name. Returns None if no attribute with the specified name exists.
- Parameters:
name (str) – The name of the attribute.
- Returns:
The attribute with the specified name or None.
- Return type:
- get_attribute_interpolator(attribute_name)
Returns interpolator for attribute on pathline, or None if it does not exist.
- Parameters:
attribute_name (str|None) – string indicating name of attribute
- Returns:
interpolator for attribute on pathline, or None if it does not exist
- Return type:
interp1d|None
- get_attribute_names()
Returns a list of all attributes stored in pathline.
- Returns:
A list of attribute names.
- Return type:
List[str]
- get_position_attribute()
Returns the ‘Position’ attribute of the pathline. Has to exist by definition.
- Raises:
AttributeError – If no position attribute exists on pathline.
- Returns:
The ‘Position’ attribute.
- Return type:
- get_reason_for_termination()
Returns the reason for termination of the pathline integration.
- Returns:
The reason for termination.
- Return type:
str
- get_t0()
Returns the initial integration time of the pathline.
- Returns:
The initial integration time.
- Return type:
float
- get_tend()
Returns the final integration time of the pathline.
- Returns:
The final integration time.
- Return type:
float
- unify_attributes(attribute_names, ref_attribute_name='Position')
Interpolates all attributes to the same time points of some reference attribute (if not specified, use the ‘Position’ attribute) and return them as an array. The first column of each array contains the time points used.
- Parameters:
attribute_names (List[str] | None) – A list of attribute names to unify. If None, all attributes are unified.
ref_attribute_name (str) – The name of the attribute to use as reference, defaults to ‘Position’.
- Returns:
A tuple containing the array of interpolated values and a list of attribute names. The first column of the array contains the time points used.
- Return type:
Tuple[NDArray, List[str]]
- class hemtracer.pathlines.PathlineAttribute(t, val, name, interpolation_scheme='linear')
Bases:
object
Class for storing a quantity of interest along a pathline. The data is stored as a function of the integration time. All class attributes are meant to be publicly accessible.
- Parameters:
t (ArrayLike) – The integration times along the pathline.
val (ArrayLike) – The values of the attribute along the pathline.
name (str) – The name of the attribute.
interpolation_scheme (str) – The scheme to use for interpolation, defaults to ‘linear’. Refer to scipy.interpolate.interp1d for other options.
- get_number_of_components()
Returns the number of components of the attribute.
- Returns:
The number of components.
- Return type:
int
- interpolator
A function that interpolates the values along the pathline.
- name
The name of the quantity of interest.
- t
The integration times along the pathline at which the values are stored.
- val
The values of the quantity of interest along the pathline.
- class hemtracer.pathlines.PathlineCollection
Bases:
object
Class for storing a collection of pathlines.
- _write_csv(filename, pathlines)
Write pathlines to CSV file. It is assumed that all pathlines contain the same attributes.
- Parameters:
filename (str) – The name of the file to write to, including path and extension.
pathlines (List[Tuple[NDArray,List[str]]]) – A list of tuples, each containing the array of interpolated values and a list of attribute names of a pathline. The first column of the array contains the time points used.
- Return type:
None
- _write_npz(filename, pathlines)
Write pathlines to npz file. It is not assumed that all pathlines contain the same attributes. The attribute names and values on each pathline are stored as separate arrays in the file, pattern: arr_0 = attribute_values_pl1, arr_1 = attribute_names_pl1, arr_2 = attribute_values_pl2, arr_3 = attribute_names_pl2, ….
- Parameters:
filename (str) – The name of the file to write to.
pathlines (List[Tuple[NDArray,List[str]]]) – A list of tuples, each containing the array of interpolated values and a list of attribute names of a pathline. The first column of the array contains the time points used.
- Return type:
None
- _write_vtk(filename, pathlines)
Write pathlines to VTK file. All pathlines are written to the same VTK file. They are identified by the attribute ‘PathlineID’, which is stored as a point data array.
- Parameters:
filename (str) – The name of the file to write to.
pathlines (List[Tuple[NDArray,List[str]]]) – A list of tuples, each containing the array of interpolated values and a list of attribute names of a pathline. The first column of the array contains the time points used.
- Return type:
None
- get_attribute(attribute_name)
Returns a list of dictionaries, each one representing a pathline and containing the keys ‘t’ and ‘y’ for time and attribute, respectively.
- Parameters:
attribute_name (str) – The name of the attribute to return.
- Returns:
A list of dictionaries.
- Return type:
List[Dict[str, NDArray]]
- get_name_distance_center()
Returns the name of the attribute that describes the orthogonal distance to the center of rotation. None if not available.
- Returns:
The name of the attribute.
- Return type:
str | None
- get_name_omega_frame()
Returns the name of the attribute that describes the angular velocity vector of the frame of reference. None if not available.
- Returns:
The name of the attribute.
- Return type:
str | None
- get_name_velocity()
Returns the name of the velocity field used for pathline integration.
- Returns:
The name of the velocity field.
- Return type:
str
- get_name_velocity_gradient()
Returns the name of the attribute that contains the velocity gradient. None if not available.
- Returns:
The name of the attribute.
- Return type:
str | None
- get_pathlines()
Returns the current list of pathlines.
- Returns:
A list of pathlines.
- Return type:
List[Pathline]
- to_file(filename, attribute_names=None)
Write pathlines to file. If attribute names are specified, only those attributes are written to file. Otherwise, all attributes are written to file.
Interpolates all attributes to the same time points (those of the ‘Position’ attribute) before writing to file.
Supports writing to CSV, npz, and VTK files. In case of CSV and VTK, it is assumed that all pathlines contain the same attributes, as this is required by the file format (every point must have the same number of attributes). If you want to write pathlines with different attributes to file, there are two options: (1) Write to npz file, which does not require all pathlines to have the same attributes. (2) Instantiate multiple PathlineTracker objects, each with a different set of attributes, and write each one to a separate file. In case of npz, the attribute names and values on each pathline are stored as separate arrays in the file, pattern: arr_0 = attribute_values_pl1, arr_1 = attribute_names_pl1, arr_2 = attribute_values_pl2, arr_3 = attribute_names_pl2, ….
- Parameters:
filename (str) – The name of the file to write to, including path and extension. Supported extensions are .csv, .npz, and .vtk.
attribute_names (List[str], optional) – A list of attribute names to write to file. If None, all attributes are written to file (default).
- Return type:
None
- class hemtracer.pathlines.PathlineReader(filename, id_name, t_name, posX_name, posY_name, posZ_name, velX_name=None, velY_name=None, velZ_name=None, dvX_dx_name=None, dvX_dy_name=None, dvX_dz_name=None, dvY_dx_name=None, dvY_dy_name=None, dvY_dz_name=None, dvZ_dx_name=None, dvZ_dy_name=None, dvZ_dz_name=None, omegaX_name=None, omegaY_name=None, omegaZ_name=None, distance_center_name=None, idx=None, gradient_interp='previous')
Bases:
PathlineCollection
Class for reading pathlines from file. Currently only supports reading from csv files.
- Parameters:
filename (str) – The name of the file to read from.
id_name (str) – The name of the attribute containing the pathline IDs.
t_name (str) – The name of the attribute containing the integration times.
posX_name (str) – The name of the attribute containing the x-coordinates of the pathlines.
posY_name (str) – The name of the attribute containing the y-coordinates of the pathlines.
posZ_name (str) – The name of the attribute containing the z-coordinates of the pathlines.
velX_name (str | None) – The name of the attribute containing the x-component of the velocity field. Defaults to None.
velY_name (str | None) – The name of the attribute containing the y-component of the velocity field. Defaults to None.
velZ_name (str | None) – The name of the attribute containing the z-component of the velocity field. Defaults to None.
dvX_dx_name (str | None) – The name of the attribute containing the derivative of v_x w.r.t. x. Defaults to None.
dvX_dy_name (str | None) – The name of the attribute containing the derivative of v_x w.r.t. y. Defaults to None.
dvX_dz_name (str | None) – The name of the attribute containing the derivative of v_x w.r.t. z. Defaults to None.
dvY_dx_name (str | None) – The name of the attribute containing the derivative of v_y w.r.t. x. Defaults to None.
dvY_dy_name (str | None) – The name of the attribute containing the derivative of v_y w.r.t. y. Defaults to None.
dvY_dz_name (str | None) – The name of the attribute containing the derivative of v_y w.r.t. z. Defaults to None.
dvZ_dx_name (str | None) – The name of the attribute containing the derivative of v_z w.r.t. x. Defaults to None.
dvZ_dy_name (str | None) – The name of the attribute containing the derivative of v_z w.r.t. y. Defaults to None.
dvZ_dz_name (str | None) – The name of the attribute containing the derivative of v_z w.r.t. z. Defaults to None.
omegaX_name (str | None) – The name of the attribute containing the x-component of the angular velocity vector of the frame of reference. Defaults to None.
omegaY_name (str | None) – The name of the attribute containing the y-component of the angular velocity vector of the frame of reference. Defaults to None.
omegaZ_name (str | None) – The name of the attribute containing the z-component of the angular velocity vector of the frame of reference. Defaults to None.
distance_center_name (str | None) – The name of the attribute containing the orthogonal distance to the center of rotation. Defaults to None.
idx (List[int] | None) – A list of indices indicating which pathlines to read from file. If None, all pathlines are read. Defaults to None.
gradient_interp (str) – The interpolation scheme to use for the velocity gradient field. Defaults to ‘previous’.
- _read_csv(filename, id_name, t_name, pos_names, vel_names, vel_grad_names, omega_names, distance_center_name, idx, gradient_interp='previous')
Read pathlines from CSV file.
- Parameters:
filename (str) – The name of the file to read from.
id_name (str) – The name of the attribute containing the pathline IDs.
t_name (str) – The name of the attribute containing the integration times.
pos_names (List[str]) – A list of attribute names containing the x, y, and z-coordinates of the pathlines.
vel_names (List[str|None]) – A list of attribute names containing the x, y, and z-components of the velocity field. If None, no velocity field is read.
vel_grad_names (List[str|None]) – A list of attribute names containing the derivatives of the velocity field. If None, no velocity gradient field is read.
omega_names (List[str|None]) – A list of attribute names containing the components of the angular velocity vector of the frame of reference. If None, no angular velocity vector is read.
distance_center_name (str | None) – The name of the attribute containing the orthogonal distance to the center of rotation. If None, no distance to center is read.
idx (List[int] | None) – A list of indices indicating which pathlines to read from file. If None, all pathlines are read.
gradient_interp (str) – The interpolation scheme to use for the velocity gradient field. Defaults to ‘previous’.
- Return type:
None
- class hemtracer.pathlines.PathlineTracker(flow_field)
Bases:
PathlineCollection
Class for generating pathlines from a given flow field. Uses VTK’s vtkStreamTracer to compute pathlines. Stores pathlines as
Pathline
objects. All point-centered data available in the Eulerian field is interpolated to the pathlines. Cell-centered data is not interpolated. Data can be interpolated to the pathlines afterwards using theinterpolate_to_pathlines()
function.Associate Eulerian flow field data and find appropriate velocity field.
- Parameters:
flow_field (EulerianFlowField) – The flow field in which to track pathlines.
- compute_pathlines(x0, initial_step=0.001, min_step=0.001, max_step=0.002, max_err=0.001, max_length=5.0, n_steps=100000, terminal_velocity=1e-10, integrator='RK45', direction='forward')
Compute pathlines starting from a list of initial points. Stores pathlines in internal list. All point-centered data available in the Eulerian field is interpolated to the pathlines. Cell-centered data is not interpolated. Data can be interpolated to the pathlines afterwards using the interpolate_to_pathlines function.
- Parameters:
x0 (List[Vector3]) – A list of seed points.
initial_step (float) – The initial step size for the pathline integration. Defaults to 0.001.
min_step (float) – The minimum step size for the pathline integration. Defaults to 0.001.
max_step (float) – The maximum step size for the pathline integration. Defaults to 0.002.
max_err (float) – The maximum error for the pathline integration. Defaults to 1e-3.
max_length (float) – The maximum length of the pathline. Defaults to 5.0.
n_steps (float) – The maximum number of steps to take. Defaults to 100000.
terminal_velocity (float) – The velocity at which to stop integration. Defaults to 1e-10.
integrator (str) – The integrator to use. Options are “RK2”, “RK4” and “RK45”. Defaults to “RK45”.
direction (str) – The direction of integration. Options are “forward”, “backward” and “both”. Defaults to “forward”.
- Return type:
None
- get_flow_field()
Returns the flow field associated with the pathline tracker.
- Returns:
The flow field.
- Return type:
- get_name_distance_center()
Returns the name of the attribute that describes the orthogonal distance to the center of rotation. Implemented by asking the flow field.
- Returns:
The name of the attribute.
- Return type:
str
- get_name_omega_frame()
Returns the name of the attribute that describes the angular velocity vector of the frame of reference. Implemented by asking the flow field.
- Returns:
The name of the attribute.
- Return type:
str
- get_name_velocity_gradient()
Returns the name of the attribute that contains the velocity gradient. Implemented by asking the flow field.
- Returns:
The name of the attribute.
- Return type:
str
- interpolate_dv_to_pathlines(sampling_rate=0.001, interpolation_scheme='previous')
Special case of interpolate_to_pathlines: interpolate velocity gradient data to pathlines. If velocity gradient has not been computed yet, it is computed as part of this function.
- Parameters:
sampling_rate (float) – The sampling rate for interpolation. Positive values are interpreted as the (fixed) time interval between subsequent samples on the pathline. Zero means to use the same time points as the pathline integration. Negative values are not allowed. Defaults to 0.001.
interpolation_scheme (str) – The interpolation scheme to use along the pathline. For all possible values, refer to scipy.interpolate.interp1d. Defaults to ‘previous’.
- Return type:
None
- interpolate_to_pathlines(field_names, sampling_rate=0.0, interpolation_scheme='linear')
Interpolate field data (cell-centered or point-centered) to pathlines. If fields of same names already exist on pathlines, they are overwritten.
- Parameters:
field_names (List[str]) – A list of field names to interpolate.
sampling_rate (float) – The sampling rate for interpolation. Positive values are interpreted as the (fixed) time interval between subsequent samples on the pathline. Zero means to use the same time points as the pathline integration. Defaults to 0.0.
interpolation_scheme (str) – The interpolation scheme to use along the pathline. For all possible values, refer to scipy.interpolate.interp1d. Defaults to ‘linear’.
- Return type:
None
hemolysis_solver
- class hemtracer.hemolysis_solver.HemolysisSolver(pathlines)
Bases:
object
Class for computing hemolysis along pathlines. Takes an existing pathline collection and handles the various pathlines contained within, interfacing to the hemolysis models.
Associate pathlines with hemolysis solver and get names of relevant quantities.
- Parameters:
pathlines (PathlineCollection) – Pathlines to analyze.
- average_hemolysis(model, end_criterion_attribute=None, end_criterion_value=0.0)
Average hemolysis index over the end points of all pathlines.
- Parameters:
model (PowerLawModel) – Power law model to use.
end_criterion_name – Criterion for selecting end points. Can be a string representing an attribute name. If given, uses the first time point where the attribute is greater than end_criterion_value as end point. Discards all pathlines where the criterion is not met. Defaults to None.
end_criterion_value (float) – Value to use as end criterion. Defaults to 0.0.
end_criterion_attribute (str | None) –
- Returns:
Average hemolysis index.
- Return type:
float
- compute_hemolysis(powerlaw_model, visc=0.0035)
Computes index of hemolysis along pathlines in percent.
- Parameters:
powerlaw_model (PowerLawModel) – Power law model to use for computing index of hemolysis.
visc (float | str) – Blood viscosity. Defaults to 0.0035 Pa s. If a string is given, it is assumed to be the name of an attribute containing the (local) dynamic viscosity.
- Return type:
None
- compute_representativeShear(model, store_solution=False)
Obtain representative scalar shear rate (effective shear) from stress-based or strain-based cell model.
- Parameters:
model (RBCModel) – Cell model to use.
store_solution (bool) – Store solution in pathlines. Defaults to False.
- Return type:
None
- get_output(model)
Obtain hemolysis solutions along pathlines after they have been computed. Returns a list of dictionaries, each one representing a pathline and containing the keys ‘t’ and ‘y’ for time and output variable.
- Parameters:
model (str) – Model to consider.
- Returns:
List of dictionaries, each one representing a pathline and containing the keys ‘t’ and ‘y’ for time and output variable.
- Return type:
List[Dict[str, NDArray]]
rbc_model
This module defines cell models used in hemolysis computations, i.e., stress-based and strain-based models for cell deformation.
- class hemtracer.rbc_model.rbc_model.RBCModel
Abstract base class for any model that computes a scalar shear rate along a pathline. These models assume a certain red blood cell (RBC) behavior, hence the name RBC Model. They can be categorized into two groups: stress-based and strain-based models. Stress-based models assume direct deformation of the cell in response to stress, reducing the instantaneous three-dimensional local fluid strain to a representative scalar shear rate. Strain-based models take into account the characteristic membrane relaxation time, employing a differential equation to explicitly resolve cell deformation in response to the acting flow forces. They then compute a scalar shear rate (so-called effective shear) from cell deformation.
- _compute_second_invariant(T)
Compute second invariant of tensor.
- Parameters:
T (Tensor3) – Tensor.
- Returns:
Second invariant.
- Return type:
float
- _compute_shear_rate_second_invariant(E)
Compute shear rate from strain tensor using second invariant.
- Parameters:
E (Tensor3) – Strain tensor.
- Returns:
Shear rate.
- Return type:
float
- _compute_strain_tensor(t)
Compute strain tensor at time t.
- Parameters:
t (float) – Time.
- Returns:
Strain tensor.
- Return type:
Tensor3
- compute_representative_shear()
Solve blood damage model to compute scalar shear rate and return times and scalar shear rate along pathline. The third return can be used to return the solution of ODE models.
- Returns:
Tuple of times and scalar shear rate along pathline.
- Return type:
Tuple[NDArray, NDArray, NDArray|None]
- get_attribute_name()
Get name of resulting scalar shear rate attribute on pathlines.
- Returns:
Name of scalar shear rate attribute.
- Return type:
str
- get_name()
Get name of blood damage model. Should not contain spaces. Has to be implemented in subclasses.
- Returns:
Name of blood damage model.
- Return type:
str
- set_time_dependent_quantitites(t0, tend, time_points=None, dv=None, omega=None, x=None, v=None, init=None)
Set time-dependent quantities for scalar shear rate model and check if all required quantities are there. Start and end time as well as velocity gradient are always required. If omega is not defined, it is set to zero, assuming a stationary frame of reference. The tank-treading model with pathline additionally requires x and v. This is handled by
hemtracer.hemolysis_solver.HemolysisSolver
.- Parameters:
t0 (float) – Start time.
tend (float) – End time.
time_points (NDArray) – Time points at which to sample velocity gradients.
dv (Callable[[float],Vector9]) – Velocity gradient tensor (in vector form) as a function of pathline time.
omega (Callable[[float],Vector3]) – MRF angular velocity vector as a function of pathline time.
x (Callable[[float],Vector3]) – Relevant position measure as a function of pathline time, e.g., absolute coordinates or distance from rotational center. Only required for tank-treading with pathline correction.
v (Callable[[float],Vector3]) – Velocity as a function of pathline time. Only required for tank-treading with pathline correction.
init (Dict[str, float]) – Dict with initial values for attributes on pathline, only required for specific shear initial condition.
- Return type:
None
- class hemtracer.rbc_model.stress_based.stress_based.Bludszuweit
Bases:
StressBasedModel
Represents the stress-based Bludszuweit model.
Initialize stress-based model. By default, use same time points as pathline, i.e., do not sample.
- _compute_representative_shear(E)
Reduce 3D strain tensor to representative scalar shear rate using the law proposed by Bludszuweit (1995)
- Parameters:
E (Tensor3) – 3D strain tensor.
- Returns:
Scalar shear rate.
- Return type:
float
- get_name()
Get name of blood damage model.
- Returns:
Name of the blood damage model.
- Return type:
str
- class hemtracer.rbc_model.stress_based.stress_based.FaghihSharp
Bases:
StressBasedModel
Represents the stress-based Faghih and Sharp model that weighs extensional and shear stresses differently.
Initialize stress-based model. By default, use same time points as pathline, i.e., do not sample.
- _compute_representative_shear(E)
Reduce 3D strain tensor to representative scalar shear rate using the law proposed by Faghih and Sharp.
- Parameters:
E (Tensor3) –
- Return type:
float
- get_name()
Get name of blood damage model.
- Returns:
Name of the blood damage model.
- Return type:
str
- class hemtracer.rbc_model.stress_based.stress_based.Frobenius
Bases:
StressBasedModel
Computes a representative scalar from instantaneous fluid strain using the Frobenius norm.
Initialize stress-based model. By default, use same time points as pathline, i.e., do not sample.
- _compute_representative_shear(E)
Reduce 3D strain tensor to representative scalar shear rate using Frobenius norm on strain tensor.
- Parameters:
E (Tensor3) – 3D strain tensor.
- Returns:
Scalar shear rate.
- Return type:
float
- get_name()
Get name of blood damage model.
- Returns:
Name of the blood damage model.
- Return type:
str
- class hemtracer.rbc_model.stress_based.stress_based.SecondInvariant
Bases:
StressBasedModel
Computes a representative scalar from instantaneous fluid strain using the Second strain invariant.
Initialize stress-based model. By default, use same time points as pathline, i.e., do not sample.
- _compute_representative_shear(E)
Reduce 3D strain tensor to representative scalar shear rate using second strain invariant.
- Parameters:
E (Tensor3) – 3D strain tensor.
- Returns:
Scalar shear rate.
- Return type:
float
- get_name()
Get name of blood damage model.
- Returns:
Name of the blood damage model.
- Return type:
str
- class hemtracer.rbc_model.stress_based.stress_based.StressBasedModel
Bases:
RBCModel
Represents a stress-based blood damage model. These models reduce the instantaneous three-dimensional local fluid strain \(\mathbf{E}\) to a representative scalar shear rate \(G_s\). They assume that RBCs deform instantly in response to changes in fluid stress. The local velocity gradient tensor is taken from pathlines or sampled at fixed time intervals and the representative scalar shear rate is computed from the strain tensor.
Initialize stress-based model. By default, use same time points as pathline, i.e., do not sample.
- _compute_representative_shear(E)
Reduce 3D strain tensor to representative scalar shear rate.
- Parameters:
E (Tensor3) – 3D strain tensor.
- Returns:
Scalar shear rate.
- Return type:
float
- compute_representative_shear()
Compute representative scalar shear rate from flow field. Called by
HemolysisSolver
.- Returns:
Tuple containing time steps and representative scalar shear rate and None, as no ODE solution is computed.
- Return type:
Tuple[NDArray, NDArray, None]
- set_sampling_rate(sampling_rate)
Set sampling rate of velocity gradients along pathline for stress-based model. If sampling rate is positive, it is interpreted as the time interval between two consecutive samples. If sampling rate is negative, it is interpreted as the number of samples per time interval.
- Parameters:
sampling_rate (float) – Sampling rate.
- Return type:
None
- class hemtracer.rbc_model.strain_based.strain_based.AroraFullEig
Bases:
MorphologyEigFormulation
Represents the eigenvalue-eigenvector formulation of the full Eulerian formulation of the Arora model (see Full Arora model)
Arora coefficients (Arora et al. 2004) are used by default.
- _RHS(t, y)
Compute the right-hand side of the ODE.
- Parameters:
t (float) – Time.
y (Vector12) – Current state.
- Returns:
Right-hand side of the ODE.
- Return type:
Vector12
- _dQdt(lamb, Et, Wt, Q)
Computes source term for orientation tensor Q.
- Parameters:
lamb (Vector3) – Eigenvalues.
Et (Tensor3) – Transformed strain tensor.
Wt (Tensor3) – Transformed vorticity tensor.
Q (Tensor3) – Eigenvectors.
- Returns:
Source term for orientation tensor.
- Return type:
Tensor3
- _initial_condition_shear(G)
Set initial condition to steady shear state.
- Parameters:
G (float) – Shear rate.
- Returns:
Steady shear state, expressed in model variables.
- Return type:
Vector12
- _initial_condition_undeformed()
Set initial condition for undeformed cells. There is some singularity in the model if the initial condition is exactly undeformed. Therefore, we add a small perturbation.
- Returns:
Initial condition for undeformed cells with a small perturbation to avoid singularity in the model.
- Return type:
Vector12
- _orthonormalize(EV)
Orthonormalizes eigenvectors.
- Parameters:
EV (Tensor3) – Eigenvectors.
- Returns:
Orthonormalized eigenvectors.
- Return type:
Tensor3
- get_name()
Get name of blood damage model.
- Returns:
Name of the blood damage model.
- Return type:
str
- class hemtracer.rbc_model.strain_based.strain_based.AroraSimplified
Bases:
MorphologyTensorFormulation
Represents the simplified Eulerian morphology model (see Simplified Arora model)
Arora coefficients (Arora et al. 2004) are used by default.
- _RHS(t, y)
Compute the right-hand side of the ODE for the cell model.
- Parameters:
t (float) – Time.
y (NDArray) – Current state.
- Returns:
Right-hand side of the ODE.
- Return type:
NDArray
- get_name()
Get name of blood damage model.
- Returns:
Name of the blood damage model.
- Return type:
str
- class hemtracer.rbc_model.strain_based.strain_based.MorphologyCoefficients(f1, f2, f3, f2t, f3t)
Bases:
tuple
Create new instance of MorphologyCoefficients(f1, f2, f3, f2t, f3t)
- _asdict()
Return a new dict which maps field names to their values.
- _field_defaults = {}
- _fields = ('f1', 'f2', 'f3', 'f2t', 'f3t')
- classmethod _make(iterable)
Make a new MorphologyCoefficients object from a sequence or iterable
- _replace(**kwds)
Return a new MorphologyCoefficients object replacing specified fields with new values
- f1
Alias for field number 0
- f2
Alias for field number 1
- f2t
Alias for field number 3
- f3
Alias for field number 2
- f3t
Alias for field number 4
- class hemtracer.rbc_model.strain_based.strain_based.MorphologyEigFormulation
Bases:
StrainBasedModel
Abstract superclass for all morphology models derived from the Arora model that employ a formulation explicit in the eigenvalues \((\lambda_1, \lambda_2, \lambda_3)\)
Arora coefficients (Arora et al. 2004) are used by default.
- _compute_Geff_from_sol(sol)
Compute effective shear rate from solution of strain-based model.
- Parameters:
sol (NDArray (n x ndf)) – Solution of strain-based model, assumed to be written as vector in form \((\lambda_1, \lambda_2, \lambda_3, \dots)\).
- Returns:
Effective shear rate.
- Return type:
NDArray (n x 1)
- _dLdt(lamb, Et)
Computes source term for eigenvalues.
- Parameters:
lamb (Vector3) – Eigenvalues.
Et (Tensor3) – Transformed strain tensor.
- Returns:
Source term for eigenvalues.
- Return type:
Vector3
- class hemtracer.rbc_model.strain_based.strain_based.MorphologyModelCoefficients(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)
Bases:
MorphologyCoefficients
,Enum
Structure that contains morphology model coefficients \((f_1, f_2, f_3, \tilde{f}_2, \tilde{f}_3)\). The coefficients are defined as follows:
\(f_1\) : Recovery coefficient.
\(f_2\) : Strain coefficient.
\(f_3\) : Vorticity coefficient.
\(\tilde{f}_2\): Coefficient for rotational part of strain tensor (only relevant for
AroraFullEig
).\(\tilde{f}_3\): Coefficient for rotational part of vorticity tensor (only relevant for
AroraFullEig
).Create new instance of MorphologyCoefficients(f1, f2, f3, f2t, f3t)
- ARORA = (5.0, 0.00042298, 0.00042298, 1.0, 1.0)
Arora et al. (2004)
- static _generate_next_value_(name, start, count, last_values)
Generate the next value when not given.
name: the name of the member start: the initial start value or None count: the number of existing members last_values: the list of values assigned
- _member_map_ = {'ARORA': MorphologyModelCoefficients.ARORA}
- _member_names_ = ['ARORA']
- _member_type_
alias of
MorphologyCoefficients
- _new_member_(f1, f2, f3, f2t, f3t)
Create new instance of MorphologyCoefficients(f1, f2, f3, f2t, f3t)
- _unhashable_values_ = []
- _use_args_ = True
- _value2member_map_ = {(5.0, 0.00042298, 0.00042298, 1.0, 1.0): MorphologyModelCoefficients.ARORA}
- _value_repr_()
Return a nicely formatted representation string
- class hemtracer.rbc_model.strain_based.strain_based.MorphologyTensorFormulation
Bases:
StrainBasedModel
Abstract base class that represents any morphology model derived from the Arora model that uses the morphology tensor
S = [S_11, S_22, S_33, S_12, S_13, S_23]
as solution variable.Arora coefficients (Arora et al. 2004) are used by default.
- _compute_Geff_from_sol(sol)
Compute effective shear rate from solution of morphology tensor model.
- Parameters:
sol (NDArray (n x 6)) – Solution of morphology tensor model, assumed to be written as vector in form [S_11, S_22, S_33, S_12, S_13, S_23].
- Returns:
Effective shear rate.
- Return type:
NDArray (n x 1)
- _compute_eig_from_morph(S)
Compute eigenvalues of morphology tensor.
- Parameters:
S (NDArray) – Morphology tensor solution, assumed to be written as vector in form [S_11, S_22, S_33, S_12, S_13, S_23].
- Returns:
Eigenvalues at nodes.
- Return type:
NDArray
- _initial_condition_shear(G)
Set initial condition for shear rate in morphology tensor models.
- Parameters:
G (float) – Shear rate.
- Returns:
Initial condition for shear.
- Return type:
NDArray
- _initial_condition_undeformed()
Set initial condition for undeformed cells in morphology tensor models, i.e., S = I.
- Returns:
Initial condition for undeformed cells.
- Return type:
NDArray
- _pack_morphology(S)
Pack values from symmetric morphology tensor into vector.
- Parameters:
S (NDArray) – Symmetric morphology tensor.
- Returns:
Vector of morphology tensor values.
- Return type:
NDArray
- _recovery_term(S)
Compute morphology recovery term for models that include morphology tensor.
- Parameters:
S (Tensor3) – Morphology tensor.
- Returns:
Morphology recovery term.
- Return type:
Tensor3
- _strain_term(E, S)
Compute strain source term for models that include morphology tensor.
- Parameters:
E (Tensor3) – Strain tensor.
S (Tensor3) – Morphology tensor.
- Returns:
Strain source term.
- Return type:
Tensor3
- _unpack_morphology(morph_vec)
Unpack values from vector into symmetric morphology tensor.
- Parameters:
morph_vec (NDArray) – Vector of morphology tensor values.
- Returns:
Symmetric morphology tensor.
- Return type:
Tensor3
- _vorticity_term(W, S)
Compute vorticity source term for models that include morphology tensor.
- Parameters:
W (Tensor3) – Vorticity tensor.
S (Tensor3) – Morphology tensor.
- Returns:
Vorticity source term.
- Return type:
Tensor3
- class hemtracer.rbc_model.strain_based.strain_based.StrainBasedModel
Bases:
RBCModel
Represents a strain-based blood damage model. These models employ a differential equation to explicitly resolve cell deformation in response to the acting flow forces. In the current implementation, this is constructed in particular for the Arora model (Arora et al. 2004) and all models derived from it, i.e., the simplified Eulerian model Pauli et al. [10], the full Eulerian morphology model (Dirkes et al. 2023) and the tank-treading model (Dirkes et al. 2023). In theory, however, also other strain-based models could be implemented as sub-classes. Since an ODE solver is required, these models offer additional configuration options for the solver and the initial condition.
Arora coefficients (Arora et al. 2004) are used by default.
- _RHS(t, y)
Right-hand side of ODE for strain-based model. Has to be implemented for each strain-based model.
- Parameters:
t (float) – Time.
y (NDArray) – Current state.
- Returns:
Right-hand side of ODE.
- Return type:
NDArray
- Raises:
NotImplementedError – If RHS is not implemented in subclasses.
- _compute_D_from_eig(lam)
Compute cell distortion D = (L-B)/(L+B) from eigenvalues of morphology tensor.
- Parameters:
lam (NDArray (n x 3)) – Nodal morphology eigenvalues lambda.
- Returns:
Nodal cell distortion values.
- Return type:
NDArray (n x 1)
- _compute_Geff_from_D(D)
Compute effective strain rate G_eff from cell distortion D.
- Parameters:
D (NDArray) – Nodal cell distortion values.
- Returns:
Solution field for effective shear rate.
- Return type:
NDArray
- _compute_Geff_from_sol(sol)
Compute effective shear rate from solution of strain-based model.
- Parameters:
sol (NDArray (n x ndf)) – Solution of strain-based model.
- Returns:
Effective shear rate.
- Return type:
NDArray (n x 1)
- Raises:
NotImplementedError – If effective shear rate computation is not implemented in subclasses.
- _compute_steady_state_shear(G)
Computes an approximate steady state for a simple shear flow of the given intensity. The equations for steady state in shear flow with f2 != f3 have been derived by myself on a piece of paper. They are not published anywhere. The equations for f2 = f3 are available in Arora et al. (2004).
- Parameters:
G (float) – Shear rate.
- Returns:
Tuple of eigenvalues of morphology tensor and basis vectors of morphology tensor.
- Return type:
Tuple[Vector3, Tensor3]
- _compute_strain_vort_rot(t)
Compute strain, vorticity, and rotation at time t.
- Parameters:
t (float) – Time.
- Returns:
Tuple of strain tensor, vorticity tensor, and angular velocity of frame of reference.
- Return type:
Tuple[Tensor3, Tensor3, Tensor3]
- _get_shear_rate_flow(t)
Compute shear rate from flow field at time t.
- Parameters:
t (float) – Time.
- Returns:
Shear rate.
- Return type:
float
- _initial_condition_shear(G)
Set initial condition for specific shear rate. Has to be implemented by subclasses to be usable.
- Parameters:
G (float) – Specific shear rate.
- Returns:
A vector of values describing a cell in the model.
- Return type:
NDArray
- Raises:
NotImplementedError – If specific shear initial condition is not implemented in this cell model.
- _initial_condition_shearRateAttribute(attribute)
Set initial condition for specific shear rate from a pathline attribute.
- Parameters:
G (str) – Name of attribute that describes shear rate to use at initial time.
attribute (str) –
- Returns:
A vector of values describing a cell in the model.
- Return type:
NDArray
- _initial_condition_steadyShear()
Set initial condition for steady shear. Has to be implemented by subclasses to be usable.
- Returns:
A vector of values describing a cell in steady shear in the model.
- Return type:
NDArray
- Raises:
NotImplementedError – If steady shear initial condition is not implemented in this cell model.
- _initial_condition_undeformed()
Set initial condition for undeformed cells. This is the default initial condition for all strain-based models. Has to be implemented by subclasses to be usable.
- Returns:
A vector of values describing an undeformed cell in the model.
- Return type:
NDArray
- Raises:
NotImplementedError – If undeformed initial condition is not implemented in this cell model.
- _unpack_antisymmetric(vec)
Unpack values of antisymmetric vector into tensor.
- Parameters:
vec (Vector3) – Vector to unpack.
- Returns:
Antisymmetric tensor.
- Return type:
Tensor3
- compute_representative_shear()
Solve strain-based model, return effective shear rate and chosen time steps. Called by
hemtracer.hemolysis_solver.HemolysisSolver
.- Returns:
Time steps and Effective shear rate and model solution.
- Return type:
Tuple[NDArray, NDArray, NDArray]
- configure_ODE_solver(method='RK45', atol=1e-06, rtol=1e-05, first_step=0.0001, max_step=0.1)
Configure ODE solver for strain-based model. For details on the parameters, see scipy.integrate.solve_ivp
- Parameters:
method (str) – Method for ODE solver. Defaults to ‘RK45’.
atol (float) – Absolute tolerance for ODE solver. Defaults to 1e-6.
rtol (float) – Relative tolerance for ODE solver. Defaults to 1e-5.
first_step (float) – Initial step size for ODE solver. Defaults to 0.0001.
max_step (float) – Maximum step size for ODE solver. Defaults to 0.1.
- Return type:
None
- set_coefficients(coeffs)
Set coefficients to use for morphology model. If not called, the Arora coefficients (Arora et al. 2004) are used by default.
- Parameters:
coeffs (MorphologyModelCoefficients) – Coefficients for morphology model, has to be a member of the
Enum
classMorphologyModelCoefficients
, e.g.,hemtracer.rbc_model.strain_based.MorphologyModelCoefficients.ARORA
.- Return type:
None
- set_initial_condition(type='undeformed', value=None)
Set initial condition for strain-based model.
- Parameters:
type (str) – Type of initial condition, i.e., shape of cells at the start of the pathline. Currently supported: undeformed, steadyShear, specificShear, solutionVector. ‘undeformed’ represents an undeformed, i.e., spherical cell. ‘steadyShear’ represents the steady state for a simple shear flow. The shear rate is computed using the second invariant of the strain tensor at the initial position. ‘specificShear’ requires another input in value that can be either a string specifying attribute to take the shear rate from, e.g., ‘Geff’, or a scalar floating point value that specifies the shear rate. It then deforms the cell to correspond to that particular scalar shear rate at the beginning. ‘solutionVector’ allows the user to specify an initial condition directly for the solution variables of the ODE. Defaults to ‘undeformed’.
val (str|float|NDArray|None) – Value for initial condition, if applicable. Defaults to None.
value (str | float | ndarray[Any, dtype[_ScalarType_co]] | None) –
- Return type:
None
- class hemtracer.rbc_model.strain_based.strain_based.TankTreading
Bases:
MorphologyEigFormulation
Represents the tank-treading cell deformation model (see Tank-treading model),
Arora coefficients (Arora et al. 2004) are used by default.
- _RHS(t, y)
Compute the right-hand side of the ODE for the tank-treading model.
- Parameters:
t (float) – Time.
y (Vector3) – Current state.
- Returns:
Right-hand side of the ODE.
- Return type:
Vector3
- _getIdxForAxis(iax)
Get indices for projection of tensor on axis iax.
- Parameters:
iax (int) – Axis.
- Returns:
Two indices (i,j) to project tensor on axis iax.
- Return type:
Tuple[int, int]
- _getRotMat(iax, theta)
Get rotation matrix for rotation around axis iax by angle theta.
- Parameters:
iax (int) – Axis.
theta (float) – Angle.
- Returns:
Rotation matrix.
- Return type:
Tensor3
- _getSteadyOrientation2D(Eii, Eij, Ejj, Wij, li, lj)
Computes the steady orientation for a 2D ellipse.
- Parameters:
Eii (float) – Diagonal element of (symmetric) strain tensor.
Eij (float) – Off-diagonal element of (symmetric) strain tensor.
Ejj (float) – Other diagonal element of (symmetric) strain tensor.
Wij (float) – Off-diagonal element of (antisymmetric) vorticity tensor.
li (float) – Eigenvalue of morphology tensor.
lj (float) – Other eigenvalue of morphology tensor.
- Returns:
Angle for steady orientation. None, if no steady orientation exists.
- Return type:
float | None
- _getSteadyOrientation3D(lamb, E, W)
Computes the steady orientation for a 3D ellipsoid. The algorithm is presented in Dirkes et al. (2023). Warning: in contrast to the paper, we employ the opposite order of eigenvalues, i.e., lambda_1 <= lambda_2 <= lambda_3. This can lead to some index differences compared to the algorithm in the paper.
- Parameters:
lamb (Vector3) – Eigenvalues of morphology tensor.
E (Tensor3) – Strain tensor.
W (Tensor3) – Vorticity tensor.
- Returns:
Orientation tensor.
- Return type:
Tensor3
- _initial_condition_shear(G)
Set initial condition to steady shear state.
- Parameters:
G (float) – Shear rate.
- Returns:
Steady shear state, expressed in model variables.
- Return type:
Vector3
- _initial_condition_undeformed()
Set initial condition for undeformed cells.
- Returns:
Initial condition for undeformed cells.
- Return type:
Vector3
- configure_steadyOrientation_solver(tol=1e-05, maxIt=100)
Set parameters for steady state orientation computation.
- Parameters:
tol (float) – Tolerance for convergence. Defaults to 1e-5.
maxIt (int) – Maximum number of iterations. Defaults to 100.
- Return type:
None
- get_name()
Get name of blood damage model.
- Returns:
Name of the blood damage model.
- Return type:
str
- class hemtracer.rbc_model.strain_based.strain_based.TankTreadingRotationCorrection
Bases:
TankTreading
Represents the tank-treading model (Dirkes et al. 2023) with a correction term for cell rotation along the pathline. This model is still experimental and not yet published.
Arora coefficients (Arora et al. 2004) are used by default.
- _compute_strain_vort_rot(t)
Overrides _compute_strain_vort_rot from superclass to include correction term for cell rotation along the pathline.
- Parameters:
t (float) – Time.
- Returns:
Strain tensor, Vorticity tensor, Angular velocity of frame of reference, including correction term.
- Return type:
Tuple[NDArray, NDArray, NDArray]
- get_name()
Get name of blood damage model.
- Returns:
Name of the blood damage model.
- Return type:
str
- set_time_dependent_quantitites(t0, tend, time_points=None, dv=None, omega=None, x=None, v=None, init=None)
Set time-dependent quantities for the cell model. Overrides superclass method to include position measure x and velocity v.
- Parameters:
t0 (float) – Start time.
tend (float) – End time.
dv (Callable[[float], Vector9]) – Time-dependent velocity gradient dv/dt.
omega (Callable[[float], Vector3]) – Time-dependent angular velocity.
x (Callable[[float], Vector3]) – Time-dependent position measure, e.g., absolute coordinates or distance from rotational center.
v (Callable[[float], Vector3]) – Time-dependent velocity.
init (Dict[str, float]) – Dict with initial values for attributes on pathline, only required for specific shear initial condition.
time_points (ndarray[Any, dtype[_ScalarType_co]] | None) –
- Raises:
ValueError – If position measure x or velocity v is not provided.
- Return type:
None
hemolysis_model
- class hemtracer.hemolysis_model.hemolysis_model.CorrelationCoefficients(A_Hb, alpha_Hb, beta_Hb)
Bases:
tuple
Named tuple for correlation coefficients. Used in
IHCorrelation
.- Parameters:
A_Hb (float) – Coefficient \(A_\mathrm{Hb}\).
alpha_Hb (float) – Coefficient \(\alpha_\mathrm{Hb}\).
beta_Hb (float) – Coefficient \(\beta_\mathrm{Hb}\).
- A_Hb
Alias for field number 0
- _asdict()
Return a new dict which maps field names to their values.
- _field_defaults = {}
- _fields = ('A_Hb', 'alpha_Hb', 'beta_Hb')
- classmethod _make(iterable)
Make a new CorrelationCoefficients object from a sequence or iterable
- _replace(**kwds)
Return a new CorrelationCoefficients object replacing specified fields with new values
- alpha_Hb
Alias for field number 1
- beta_Hb
Alias for field number 2
- class hemtracer.hemolysis_model.hemolysis_model.IHCorrelation(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)
Bases:
CorrelationCoefficients
,Enum
Represents sets of empirical parameters \((A_\mathrm{Hb}, \alpha_\mathrm{Hb}, \beta_\mathrm{Hb})\) for the correlation \(IH = A_\mathrm{Hb} \sigma^{\alpha_\mathrm{Hb}} t^{\beta_\mathrm{Hb}}\). Magnitude of \(A_{Hb}\) is chosen such that \(IH\) is in percent.
Create new instance of CorrelationCoefficients(A_Hb, alpha_Hb, beta_Hb)
- DING_BOVINE = (9.772e-05, 1.4445, 0.2076)
Ding et al. [3]
- DING_HUMAN = (3.458e-06, 2.0639, 0.2777)
Ding et al. [3]
- DING_PORCINE = (0.0006701, 1.0981, 0.2778)
Ding et al. [3]
- GESENHUES_OPTIMIZED = (0.00023212, 1.4949, 0.33)
Gesenhues et al. [6]
- GIERSIEPEN = (3.62e-05, 2.416, 0.785)
Giersiepen et al. [7]
- SONG = (1.8e-06, 1.991, 0.765)
Song et al. [11]
- ZHANG = (1.228e-05, 1.9918, 0.6606)
Zhang et al. [13]
- static _generate_next_value_(name, start, count, last_values)
Generate the next value when not given.
name: the name of the member start: the initial start value or None count: the number of existing members last_values: the list of values assigned
- _member_map_ = {'DING_BOVINE': IHCorrelation.DING_BOVINE, 'DING_HUMAN': IHCorrelation.DING_HUMAN, 'DING_PORCINE': IHCorrelation.DING_PORCINE, 'GESENHUES_OPTIMIZED': IHCorrelation.GESENHUES_OPTIMIZED, 'GIERSIEPEN': IHCorrelation.GIERSIEPEN, 'SONG': IHCorrelation.SONG, 'ZHANG': IHCorrelation.ZHANG}
- _member_names_ = ['GIERSIEPEN', 'SONG', 'ZHANG', 'DING_HUMAN', 'DING_PORCINE', 'DING_BOVINE', 'GESENHUES_OPTIMIZED']
- _member_type_
alias of
CorrelationCoefficients
- _new_member_(A_Hb, alpha_Hb, beta_Hb)
Create new instance of CorrelationCoefficients(A_Hb, alpha_Hb, beta_Hb)
- _unhashable_values_ = []
- _use_args_ = True
- _value2member_map_ = {(1.8e-06, 1.991, 0.765): IHCorrelation.SONG, (3.458e-06, 2.0639, 0.2777): IHCorrelation.DING_HUMAN, (1.228e-05, 1.9918, 0.6606): IHCorrelation.ZHANG, (3.62e-05, 2.416, 0.785): IHCorrelation.GIERSIEPEN, (9.772e-05, 1.4445, 0.2076): IHCorrelation.DING_BOVINE, (0.00023212, 1.4949, 0.33): IHCorrelation.GESENHUES_OPTIMIZED, (0.0006701, 1.0981, 0.2778): IHCorrelation.DING_PORCINE}
- _value_repr_()
Return a nicely formatted representation string
- class hemtracer.hemolysis_model.hemolysis_model.PowerLawModel(scalar_shear, hemolysis_correlation=IHCorrelation.GIERSIEPEN, mu=0.0035, integration_scheme='basic')
Bases:
object
A power law hemolysis model is a model that, given a scalar measure for fluid stress (in our case a representative shear rate), predicts the hemolysis index along a pathline. This is done by integrating an experimental correlation for the hemolysis index along the pathline. For details, see Hemolysis models.
Initialization defines all parameters to use. They cannot be changed afterwards.
- Parameters:
scalar_shear (RBCModel | str) – Model to compute scalar shear rate, or name of attribute containing representative shear rate, e.g., ‘Geff’ if available from an Eulerian simulation.
hemolysis_correlation (HemolysisCorrelation) – Hemolysis correlation to use.
mu (float | str) – Dynamic viscosity of blood. Defaults to 3.5e-3. If a string is given, it is assumed to be the name of an attribute containing the (local) dynamic viscosity.
integration_scheme (str) – Integration scheme for hemolysis correlation. Integration schemes defined according to Taskin et al. [12]. Valid options are ‘basic’ (HI1), ‘timeDiff’ (HI2), ‘linearized’ (HI3), ‘mechDose’ (HI4), ‘effTime’ (HI5). Defaults to ‘basic’.
- _compute_HI1(t, tau)
Computes index of hemolysis by directly integrating experimental correlation.
- Parameters:
t (NDArray) – Time steps.
tau (NDArray) – Scalar representative stress.
- Returns:
Hemolysis index.
- Return type:
NDArray
- _compute_HI2(t, tau)
Computes index of hemolysis by incorporating time derivative in experimental correlation.
- Parameters:
t (NDArray) – Time steps.
tau (NDArray) – Scalar representative stress.
- Returns:
Hemolysis index.
- Return type:
NDArray
- _compute_HI3(t, tau)
Computes index of hemolysis by summing linearized damage.
- Parameters:
t (NDArray) – Time steps.
tau (NDArray) – Scalar representative stress.
- Returns:
Hemolysis index.
- Return type:
NDArray
- _compute_HI4(t, tau)
Computes index of hemolysis by accumulating mechanical dose (Grigioni et al. [9]).
- Parameters:
t (NDArray) – Time steps.
tau (NDArray) – Scalar representative stress.
- Returns:
Hemolysis index.
- Return type:
NDArray
- _compute_HI5(t, tau)
Computes index of hemolysis by using virtual time step approach (Goubergrits and Affeld [8]).
- Parameters:
t (NDArray) – Time steps.
tau (NDArray) – Scalar representative stress.
- Returns:
Hemolysis index.
- Return type:
NDArray
- _integration_scheme_name
Obtain coefficients for definitions in accordance with Taskin et al. Coefficient C for hemolysis correlation \(HI = C t^\alpha \sigma^\beta\). This corresponds to the definition in Taskin et al. [12]. The internal values of this class are thus different from the definitions in
HemolysisCorrelation
, i.e., \(C = A_\mathrm{Hb}, \, \alpha = \beta_\mathrm{Hb}, \, \beta = \alpha_\mathrm{Hb}\). This is to allow for direct comparison of the numerical schemes with the formulations in Taskin et al. [12].
- _mu
Define integration scheme.
- compute_hemolysis(t, G, mu)
Compute hemolysis along pathline. Called by
HemolysisSolver
.- Parameters:
t (NDArray) – Time steps.
G (NDArray) – Scalar shear rate.
mu (NDArray) – Dynamic viscosity.
- Returns:
Hemolysis index.
- Return type:
NDArray
- get_attribute_name()
Get the name of the attribute that will be added to pathlines.
- Returns:
The name of the attribute that will be added to pathlines.
- Return type:
str
- get_name()
Get the name of the power law hemolysis model.
- Returns:
The name of the power law hemolysis model.
- Return type:
str
- get_scalar_shear_name()
Get the name of the scalar shear rate attribute.
- Returns:
The name of the scalar shear rate attribute.
- Return type:
str