PYME.Analysis.points.spherical_harmonics module¶
Estimate spherical harmonics from a point data set Initial fitting/conversions ripped 100% from David Baddeley / scipy
- class PYME.Analysis.points.spherical_harmonics.SHShell(centre=(0, 0, 0), principle_axes=((1, 0, 0), (0, 1, 0), (0, 0, 1)), axis_scaling=(1.0, 1.0, 1.0), modes=((0, 0),), coefficients=(1,))¶
Bases:
ScaledShell
Initial work on a replacement interface for ScaledShell
Goals:
can be constructed directly as well as / instead of being fit
fitting is a single function call (rather than 3)
easily serialised (TODO)
does not own/keep a copy of data points (to facilitate serialisation)
when directly constructed, can used to represent synthetic nuclei etc … whilst retaining distance evaluation functions
- fit(points, n_max=3, n_iters=2, tol=0.3, principle_axis_sampling=1.0)¶
- property scaling_factors¶
- class PYME.Analysis.points.spherical_harmonics.ScaledShell(sampling_fraction=1.0)¶
Bases:
object
- approximate_image_bounds(d_zenith=0.1, d_azimuth=0.1)¶
- approximate_normal(x, y, z, d_azimuth=1e-06, d_zenith=1e-06, return_orthogonal_vectors=False)¶
Numerically approximate a vector(s) normal to the spherical harmonic shell at the query point(s).
For input point(s), scale and convert to spherical coordinates, shift by +/- d_azimuth and d_zenith to get ‘phantom’ points in the plane tangent to the spherical harmonic expansion on either side of the query point(s). Scale back, convert to cartesian, make vectors from the phantom points (which are by definition not parallel) and cross them to get a vector perpindicular to the plane.
- Parameters
- xndarray, float
cartesian x location of point(s) on the surface to calculate the normal at
- yndarray, float
cartesian y location of point(s) on the surface to calculate the normal at
- zndarray, float
cartesian z location of point(s) on the surface to calculate the normal at
- d_azimuthfloat
azimuth step size for generating vector in plane of the shell [radians]
- d_zenithfloat
zenith step size for generating vector in plane of the shell [radians]
- Returns
- normal_vectorndarray
cartesian unit vector(s) normal to query point(s). size (len(x), 3)
- orth0ndarray
cartesian unit vector(s) in the plane of the spherical harmonic shell at the query point(s), and perpendicular to normal_vector
- orth1ndarray
cartesian unit vector(s) orthogonal to normal_vector and orth0
- check_inside(x=None, y=None, z=None)¶
- data_type = [('modes', '<2i4'), ('coefficients', '<f4')]¶
- distance_to_shell(query_points, d_angles=0.1)¶
- Parameters
- query_pointslist-like of ndarrays
Arrays of positions to query (in cartesian coordinates), i.e. [np.array(x), np.array(y), np.array(z)]
- d_anglesfloat
Sets the step size in radians of zenith and azimuth arrays used in reconstructing the spherical harmonic shell
- Returns
- min_distancefloat
minimum distance from query points (i.e. input coordinate) to the spherical harmonic surface
- closest_points_on_surfacetuple of floats
returns the position in cartesian coordinates of the point on the surface closest to the input ‘position’
- distance_to_shell_along_vector_from_point(vector, starting_point, guess=None)¶
Calculate the distance to the shell along a given direction, from a given point.
- Parameters
- vectorlist-like
cartesian vector indicating direction to query for proximity to shell
- starting_pointlist-like
cartesian position from which to start traveling along input vector when calculating shell proximity
- guess_distancesarray, float
initial guess for distance solver. See self._distance_error()
- Returns
- fit_shell(max_n_mode=3, max_iterations=2, tol_init=0.3)¶
- static from_tabular(shell_table, index=0)¶
- get_fitted_shell(azimuth, zenith)¶
- get_mesh_vertices_faces(d_zenith=0.1)¶
Compute vertices and faces to pass to PYME.experimental._triangle_mesh.TriangleMesh.
Note this will be non-manifold because of cuts at zenith=0, azimuth=0.
- Parameters
- d_zenithfloat, optional
zenith step size for generating vector in plane of the shell [radians], by default 0.1
- radial_distance_to_shell(points)¶
finds distance to shell along a vector from the centre of the shell to a point.
Notes
This is not a geometric/euclidean distance, but should be a good approximation close to the shell (or far away). Over mid-range distances, it will give potentially very distorted distances. It is, however, very much faster and can safely be used in applications where strict scaling is not important - e.g. for a SDF representation of the surface.
- set_fitting_points(x, y, z)¶
- shell_coordinates(points)¶
Scale query points, projecting them onto the basis used in shell- fitting. Return in scaled spherical coordinates
- Parameters
- points3-tuple
x, y, z positions
- Returns
- azimuth, zenith, r
scaled spherical coordinates of input points
- property table_representation¶
- to_hdf(*args, **kwargs)¶
- to_recarray(keys=None)¶
Pretend we are a PYME.IO.tabular type
- Parameters
- keysNone
Ignored for this contrived function
- Returns
- numpy recarray version of self
- uniform_random_radial_density(n_radial_bins=20, batch_size=100000, target_sampling_nm=75.0)¶
Estimate a radial density histogram for uniform spatial sampling within the shell. Performs calculation in batches to give bounded memory usage.
- Parameters
- n_radial_binsint
Number of radial bins to divide the radius into
- batch_sizeint
number of test points to generate in one batch. Note that the default (100000) might be a bit conservative
- target_sampling_nm: float
desired sampling density of combined batches. Used in conjuction with batch_size to determine the number of batches required.
- PYME.Analysis.points.spherical_harmonics.generate_icosahedron()¶
Generate an icosahedron in spherical coordinates.
- Returns
- azimuthnp.array
Length 12, azimuth of icosahedron vertices [0, 2*pi].
- zenithnp.array
Length 12, zenith of icosahedron vertices [0, pi].
- facesnp.array
20 x 3, counterclockwise triangles denoting icosahedron faces.
- PYME.Analysis.points.spherical_harmonics.icosahedron_mesh(n_subdivision=1)¶
Return an icosahedron subdivided n_subdivision times. This provides a quasi-regular sampling of the unit sphere.
v0 v0 / \ / / \ ==> v3----v5 / \ / \ / v1------v2 v1--v4--v2
- Returns
- azimuthnp.array
Azimuth of mesh vertices [0, 2*pi].
- zenithnp.array
Zenith of mesh vertices [0, pi].
- facesnp.array
20 x 3, counterclockwise triangles denoting mesh faces.
- PYME.Analysis.points.spherical_harmonics.r_sph_harm(m, n, azimuth, zenith)¶
return real valued spherical harmonics. Uses the convention that m > 0 corresponds to the cosine terms, m < zero the sine terms
- Parameters
- mint
order of the spherical harmonic, |m| <= n
- nint
degree of the spherical harmonic, n >= 0
- azimuthndarray
the azimuth angle in [0, 2pi]
- zenithndarray
the elevation in [0, pi]
- Returns
- PYME.Analysis.points.spherical_harmonics.reconstruct_shell(modes, coeffs, azimuth, zenith)¶
- PYME.Analysis.points.spherical_harmonics.scaled_shell_from_hdf(hdf_file, table_name='harmonic_shell')¶
- Parameters
- hdf_filestr or tables.file.File
path to hdf file or opened hdf file
- table_namestr
name of the table containing the spherical harmonic expansion information
- Returns
- shellScaledShell
see nucleus.spherical_harmonics.shell_tools.ScaledShell
- PYME.Analysis.points.spherical_harmonics.sphere_expansion(x, y, z, n_max=3)¶
Project coordinates onto spherical harmonics
- Parameters
- xndarray
x coordinates
- yndarray
y coordinates
- zndarray
z coordinates
- n_maxint
Maximum order to calculate to
- Returns
- modeslist of tuples
a list of the (m, n) modes projected onto
- cndarray
the mode coefficients
- PYME.Analysis.points.spherical_harmonics.sphere_expansion_clean(x, y, z, n_max=3, max_iters=2, tol_init=0.3)¶
Project coordinates onto spherical harmonics
- Parameters
- xndarray
x coordinates
- yndarray
y coordinates
- zndarray
z coordinates
- n_maxint
Maximum order to calculate to
- max_iters: int
number of fit iterations
- tol_init: float
relative outlier tolerance. Used to ignore outliers in subsequent iterations
- Returns
- modeslist of tuples
a list of the (m, n) modes projected onto
- cndarray
the mode coefficients