Measurements module#
Measurements subpackage for napari-stress.
This subpackage contains functions for measuring curvature, stress, and temporal correlations on surfaces.
- napari_stress.measurements.anisotropic_stress(mean_curvature_droplet: ndarray, H0_droplet: float, mean_curvature_ellipsoid: ndarray, H0_ellipsoid: float, gamma: float = 26.0) Tuple[ndarray, ndarray, ndarray] #
Calculate anisotropic stress from mean and averaged curvatures.
- Parameters:
mean_curvature_droplet (np.ndarray) – mean curvature at every point on the surface of a droplet
H0_droplet (float) – surface-integrated surface curvature on droplet
mean_curvature_ellipsoid (np.ndarray) – mean curvature at every point on the surface of an ellipsoid that was fitted to a droplet. The droplet positions must correspond to the point locations on the droplet surface in terms of latitude and longitude
H0_ellipsoid (float) – surface-integrated surface curvature on ellipsoid
gamma (float, optional) – interfacial surface tension in mN/m. The default is 26.0. See also [1].
- Returns:
stress (np.ndarray) – raw anisotropic stress on every point on the droplet surface
stress_tissue (np.ndarray) – tissue-scale anisotropic stress on the droplet surface
stress_droplet (np.ndarray) – cell-scale anisotropic stress on the droplet surface
References
- napari_stress.measurements.average_mean_curvatures_on_manifold(input_manifold: manifold) Tuple[float, float] #
Calculate averaged mean curvatures on manifold.
The input can also be a layer with a manifold in its metadata. In this case the results are stored in the layer’s metadata.
- Parameters:
manifold (mnfd.manifold)
- Returns:
H0_arithmetic (float) – Arithmetic average of mean curvature
H0_surface_integral (float) – Averaged curvature by integrating surface area
H0_volume_integral (float) – Averaged curvature by deriving volume-integral of manifold.
S2_volume (float) – Volume of the unit sphere
H0_radial_surface (float) – Averaged curvature on radially expanded surface. Only calculated if mnfd.manifold.manifold_type is radial.
- napari_stress.measurements.calculate_anisotropy(df: DataFrame, column: str, alpha: float = 0.05, group_column: str = 'time') DataFrame #
Calculate anisotropy of a column in a dataframe.
The dataframe is assumed to contain multiple groups, which are defined by the values in the group_column.
- Parameters:
df (pd.DataFrame) – Dataframe containing the data to analyze
column (str, optional) – Column name to analyze. The column is assumed to contain numerical data.
alpha (float, optional) – Lower and upper percentile of the data to exclude when calculating the anisotropy, by default 0.05
group_column (str, optional) – Column name to use for grouping the data, by default ‘time’
- Returns:
Dataframe containing the anisotropy of the data in the column for every group in the dataframe: - column + ‘_lower’: lower percentile of the data - column + ‘_upper’: upper percentile of the data - column + ‘_anisotropy’: anisotropy of the data
- Return type:
pd.DataFrame
- napari_stress.measurements.calculate_mean_curvature_on_manifold(input_manifold: manifold) Tuple[ndarray, float, float] #
Calculate mean curvatures for a given manifold.
- Parameters:
manifold (mnfd.manifold) – Input manifold to calculate mean curvature for.
- Returns:
mean_curvatures – Mean curvature value for every quadrature point
- Return type:
np.ndarray
- napari_stress.measurements.calculate_patch_fitted_curvature_on_pointcloud(points: napari.types.PointsData, search_radius: float = 2) DataFrame #
Calculate the curvature of a patch fitted to the points.
- Parameters:
points ('napari.types.PointsData') – The points to calculate the curvature for.
- Returns:
Dataframe with columns “mean_curvature”, “principal_curvature_1” and “principal_curvature_2” for the mean curvature and principal curvatures of the patch fitted to every point in the pointcloud.
- Return type:
pd.DataFrame
- napari_stress.measurements.calculate_patch_fitted_curvature_on_surface(surface: napari.types.SurfaceData, search_radius: float = 2) DataFrame #
Calculate the curvature of a patch fitted to the surface.
- Parameters:
surface ('napari.types.SurfaceData') – The surface to calculate the curvature for.
- Returns:
Dataframe with columns “mean_curvature”, “principal_curvature_1” and “principal_curvature_2” for the mean curvature and principal curvatures of the patch fitted to every vertex on the surface.
- Return type:
pd.DataFrame
- napari_stress.measurements.comprehensive_analysis(pointcloud: PointsData, max_degree: int = 5, n_quadrature_points: int = 110, maximal_distance: int = None, gamma: float = 26.0, verbose=False) List[LayerDataTuple] #
Comprehensive stress analysis of droplet points layer.
- Parameters:
pointcloud (PointsData) – Points layer.
max_degree (int, optional) – Maximum degree of spherical harmonics expansion, by default 5
n_quadrature_points (int, optional) – Number of quadrature points, by default 110
maximal_distance (int, optional) – Maximal distance for geodesic distance matrix, by default None
gamma (float, optional) – Surface tension, by default 26.0
verbose (bool, optional) – Show progress bar, by default False
- Returns:
List of layer data tuples: - layer_spherical_harmonics: ‘napari.types.PointsData’
fitted spherical harmonics expansion
- layer_fitted_ellipsoid_points: ‘napari.types.PointsData’
fitted ellipsoid points
- layer_fitted_ellipsoid: ‘napari.types.VectorsData’
fitted ellipsoid major axes
- layer_quadrature_ellipsoid: ‘napari.types.PointsData’
Lebedev quadrature points on ellipsoid
- layer_quadrature: ‘napari.types.PointsData’
Lebedev quadrature points on spherical-harmonics fitted droplet
- max_min_geodesics_total: ‘napari.types.SurfaceData’
geodesics on total stress from local maxima of total stress to local minima of total stress
- min_max_geodesics_total: ‘napari.types.SurfaceData’
geodesics on total stress from local minima of total stress to local maxima of total stress
- max_min_geodesics_cell: ‘napari.types.SurfaceData’
geodesics on cell stress from local maxima of cell stress to local minima of cell stress
- min_max_geodesics_cell: ‘napari.types.SurfaceData’
geodesics on cell stress from local minima of cell stress to local maxima of cell stress
- layer_surface_autocorrelation: ‘napari.types.SurfaceData’
Surface representation of autocorrelations of total stress.
- Return type:
List[LayerDataTuple]
- napari_stress.measurements.correlation_on_surface(surface1: napari.types.SurfaceData, surface2: napari.types.SurfaceData, distance_matrix: ndarray = None, maximal_distance: float = None) dict #
Calculate (auto-) correlation of features on surface.
This calculates the correlation of features on a surface with itself (auto-correlation) or with another surface (cross-correlation). If the two input surfaces are identical, this is the auto-correlation.
- Parameters:
surface1 ('napari.types.SurfaceData')
surface2 ('napari.types.SurfaceData')
- Returns:
- Dictionary with keys:
auto_correlations_distances
auto_correlations_average
auto_correlations_normalized
auto_correlations_averaged_normalized
- Return type:
dict
- napari_stress.measurements.curvature_on_ellipsoid(ellipsoid: napari.types.VectorsData, sample_points: napari.types.PointsData) napari.types.LayerDataTuple #
Calculate mean curvature on ellipsoid.
- Parameters:
ellipsoid (VectorsData) – Ellipsoid major axes to calculate mean curvature for.
sample_points (PointsData) – Sample points to calculate mean curvature on.
- Returns:
The sample points, properties and layer type. The properties contain the mean curvature, principal curvatures and the averaged mean curvature on the ellipsoid.
- Return type:
LayerDataTuple (tuple)
See also
- Mean curvature.
Mean curvature on Wikipedia.
- Ellipsoid definition
Ellipsoid definition on Wolfram MathWorld.
- napari_stress.measurements.deviation_from_ellipsoidal_mode(points: PointsData, max_degree: int = 5, viewer: Viewer = None) LayerDataTuple #
- napari_stress.measurements.distance_to_k_nearest_neighbors(points: napari.types.PointsData, k=5) DataFrame #
Calculate the distance to the k nearest neighbors for each point.
- Parameters:
points ('napari.types.PointsData') – The points to calculate the distance to the k nearest neighbors for.
k (int, optional) – The number of nearest neighbors to use for the calculation, by default 5.
- Returns:
A dataframe with the distance to the k nearest neighbors for each point.
- Return type:
pd.DataFrame
- napari_stress.measurements.gauss_bonnet_test(input_manifold: manifold, viewer: napari.Viewer = None) Tuple[float, float] #
Use Gauss-Bonnet theorem to measure resolution on manifold.
- Parameters:
input_manifold (mnfd.manifold) – Manifold to calculate Gauss-Bonnet error for.
- Returns:
Gauss_Bonnet_Err (float) – Absolute error of Gauss-Bonnet-test
Gauss_Bonnet_Rel_Err (float) – Relative error of Gauss-Bonnet test (absolute error divided by 4*pi)
See also
- napari_stress.measurements.geodesic_distance_matrix(surface: SurfaceData, show_progress: bool = False) ndarray #
Calculate a pairwise distance matrix for vertices of a surface.
This calculates the geodesic distances between any two vertices on a given surface.
- Parameters:
surface (SurfaceData)
- Returns:
distance_matrix – Triangular matrix with differences between vertex i and j located at distance_matrix[i, j]
- Return type:
np.ndarray
- napari_stress.measurements.geodesic_path(surface: SurfaceData, index_1: int, index_2: int) napari.types.VectorsData #
Calculate the geodesic path between two index-defined surface vertices .
- Parameters:
surface ('napari.types.SurfaceData')
index_1 (int) – Index of start vertex
index_2 (int) – Index of destination vertex
- Returns:
VectorsData array with shape (n_points, 2, 3) where n_points is the number of points on the path. The first dimension is the point coordinates, the second dimension is the vector from point to point.
- Return type:
napari.types.VectorsData
- napari_stress.measurements.haversine_distances(degree_lebedev: int, n_lebedev_points: int) ndarray #
Calculate geodesic (Great Circle) distance matrix from haversine formula.
- Parameters:
deg_lbdv (int) – degree of Lebedev quadrature
num_lbdv_pts (int) – number of Lebedev quadrature points
- Returns:
distance matrix of shape (num_lbdv_pts, num_lbdv_pts)
- Return type:
np.ndarray
See also
Haversine formula <https://en.wikipedia.org/wiki/Haversine_formula>
- napari_stress.measurements.local_extrema_analysis(surface: napari.types.SurfaceData, distance_matrix: ndarray = None) List[LayerDataTuple] #
Get local maximum and minimum and analyze their mutual distances.
This function identifies local maxima and minima on a given surface and analyzes the distances and value differences between them.
- Parameters:
surface (napari.types.SurfaceData) – The surface data on which local extrema are to be found.
distance_matrix (np.ndarray, optional) – Geodesic distance matrix for the surface. If not provided, it will be computed. The default is None.
- Returns:
A list of LayerDataTuples, each containing features and metadata related to local extrema. The ‘features’ key includes:
local_max_and_min: -1 if local minimum, +1 if local maximum, 0 otherwise.
- The ‘metadata’ key includes:
nearest_min_max_dists: Distance between each local maximum and its nearest local minimum.
nearest_min_max_anisotropies: Difference in input value (vertices, faces, values) between each local maximum and its nearest local minimum.
min_max_pair_distances: Distances between all pairs of local minima and maxima.
min_max_pair_anisotropies: Difference in input value (vertices, faces, values) between all pairs of local minima and maxima.
- Return type:
List[LayerDataTuple]
- napari_stress.measurements.maximal_tissue_anisotropy(ellipsoid: napari.types.VectorsData, gamma: float) float #
Calculate maximaum stress anisotropy on ellipsoid.
- Parameters:
ellipsoid ('napari.types.VectorsData')
gamma (float, optional) – Interfacial surface tension in mN/m.
- Returns:
maximal_tissue_anisotropy – Maximum stress anisotropy on ellipsoid
- Return type:
float
- napari_stress.measurements.mean_curvature_differences_radial_cartesian_manifolds(manifold_cartesian: manifold, manifold_radial: manifold) ndarray #
Calculate difference of radial and cartesian mean curvature calculation
- Parameters:
manifold_cartesian (manifold) – Cartesian manifold
manifold_radial (manifold) – Radial manifold
- Returns:
difference of mean curvatures.
- Return type:
np.ndarray
- napari_stress.measurements.mean_curvature_on_ellipse_cardinal_points(ellipsoid: napari.types.VectorsData) list #
Calculate mean points on major axes tip points of ellipsoid.
- Parameters:
ellipsoid ('napari.types.VectorsData') – Ellipsoid major axes to calculate mean curvature for.
- Returns:
Mean curvature at cardinal points on ellipsoid, e.g., at the inter- section of the ellipsoid major axes and the allipsoid surface.
- Return type:
list
- napari_stress.measurements.mean_curvature_on_radial_manifold(input_manifold: manifold) ndarray #
Calculate mean curvature on radial manifold
- Parameters:
input_manifold (manifold) – radial manifold
- Returns:
mean curvature
- Return type:
np.ndarray
- napari_stress.measurements.measure_intensity_on_surface(surface: napari.types.SurfaceData, image: napari.types.ImageData, measurement_range: float = 1.0, sampling_distance: float = 0.5, center: bool = True, interpolation_method: str = 'linear') DataFrame #
Measure intensity on surface.
- Parameters:
surface ('napari.types.SurfaceData')
image ('napari.types.ImageData')
measurement_range (float, optional) – Range of measurement, by default 1.0. This determines in which range around the surface the intensity is measured.
sampling_distance (float, optional) – Distance between samples, by default 0.5. This determines how many samples are taken in the measurement range. Needs to be smaller than measurement_range.
center (bool, optional) – Center normal vectors on surface, by default True. If set to False, the normal vectors point away from the surface so that the intensity is measured only on one side of the surface.
interpolation_method (str, optional) – Interpolation method, by default “linear”
- Returns:
Intensity values for each point on the surface.
- Return type:
pd.DataFrame
- napari_stress.measurements.radial_surface_averaged_mean_curvature(input_manifold: manifold) float #
Eestimate H0 on Radial Manifold.
Integratal of H_Rad (on surface) divided by Rad surface area
- Parameters:
input_manifold (manifold) – radial manifold
- Returns:
averaged mean curvature
- Return type:
float
- napari_stress.measurements.sample_intensity_along_vector(sample_vectors: napari.types.VectorsData, image: napari.types.ImageData, sampling_distance: float = 1.0, interpolation_method: str = 'linear') DataFrame #
Sample intensity along sample_vectors of equal length.
- Parameters:
sample_vectors ('napari.types.VectorsData') – Vectors along which to measure intensity. Must be of equal length.
image ('napari.types.ImageData') – Image to sample intensity from.
sampling_distance (float, optional) – Distance between samples, by default 1.0
interpolation_method (str, optional) – Interpolation method, by default “linear”. See scipy.interpolate.RegularGridInterpolator for available methods.
- Returns:
Intensity values for each point on the sample_vectors.
- Return type:
pd.DataFrame
See also
- napari_stress.measurements.spatio_temporal_autocorrelation(surfaces: SurfaceData, distance_matrix: ndarray, maximal_distance: float | None = None)#
Spatio-temporal autocorrelation.
- Parameters:
surfaces (SurfaceData) – 4D-surface object with values
distance_matrix (np.ndarray) – Distance matrix to be used on the surface.
maximal_distance (float, optional) – _description_. Defaults to None.
- napari_stress.measurements.temporal_autocorrelation(df: DataFrame, feature: str, frame_column_name: str = 'frame')#
Calculate temporal autocorrelation for a list of features.
- Parameters:
features (list) – List of features - each entry corresponds to features
timeframe (of a single)
- Returns:
temporal autocorrelation. The i-th entry denotes the correlation of features at time i with the feature at time 0.
- Return type:
np.ndarray
- napari_stress.measurements.tissue_stress_tensor(ellipsoid: napari.types.VectorsData, H0_ellipsoid: float, gamma: float) Tuple[ndarray, ndarray] #
Calculate tissue stress tensor(s).
- Parameters:
ellipsoid ('napari.types.VectorsData') – Ellipsoid that was fitted to a droplet.
H0_ellipsoid (float) – averaged mean curvature of the ellipsoid.
gamma (float) – droplet interfacial tension in mN/m
- Returns:
Tissue_Stress_Tensor_elliptical (np.ndarray) – 3x3 orientation matrix with stresses along ellipsoid axes
Tissue_Stress_Tensor_cartesian (np.ndarray) – 3x3 orientation matrix with stresses projected onto cartesian axes