SpectralCloudstering¶
-
class
scimes.
SpectralCloudstering
(dendrogram, catalog, header, criteria=['volume', 'luminosity'], user_k=None, user_ams=None, user_scalpars=None, user_iter=None, save_isol_leaves=False, save_clust_leaves=False, save_branches=False, save_all_leaves=False, save_all=False, blind=False, rms=<Mock name='mock.nan' id='140089526284816'>, s2nlim=3, locscaling=False)[source] [edit on github]¶ Bases:
object
Apply the spectral clustering to find the best cloud segmentation out from a dendrogram.
Parameters: - dendrogram: ‘astrodendro.dendrogram.Dendrogram’ instance
The dendrogram to clusterize.
- catalog: ‘astropy.table.table.Table’ instance
A catalog containing all properties of the dendrogram structures. Generally generated with ppv_catalog module.
- header: ‘astropy.io.fits.header.Header’ instance
The header of the fits data the dendrogram was generated from. Necessary to obtain the assignment cubes.
- criteria: list of strings
Clustering criteria referred to the structure properties in the catalog (default [‘volume’, ‘luminosity’]).
- user_k: int
The expected number of clusters, if not provided it will be guessed automatically through the eigenvalues of the unsmoothed affinity matrix.
- user_ams: numpy array
User provided affinity matrix. Whether this is not furnish it is automatically generated through the volume and/or luminosity criteria.
- user_scalpars: list of floats
User-provided scaling parameters to smooth the affinity matrices.
- user_iter: int
User-provided number of k-means iterations.
- save_isol_leaves: bool
Consider the isolated leaves (without parent) as individual ‘clusters’. Useful for low resolution data where the beam size corresponds to the size of a Giant Molecular Cloud.
- save_clust_leaves: bool
Consider unclustered leaves as individual ‘clusters’. This keyword will not include the isolated leaves without parents.
- save_all_leaves: bool
Trigger both save_isol_leaves and save_clust_leaves.
- save_branches: bool
Retain all isolated branches usually discarded by the cluster analysis.
- save_all: bool
Trigger all save_isol_leaves, save_clust_leaves, and save_branches.
- rms: int or float
Noise level of the observation. Necessary to calculate the scaling parameter above a certain signal-to-noise ratio.
- s2nlim: int or float
Signal-to-noise limit above which the scaling parameter is calculated. Needed only if rms is not np.nan.
- blind: bool
Show the affinity matrices. Matplotlib required.
- locscaling: bool
Smooth the affinity matrices using a local scaling technique. This does not work well …
Methods Summary
plot_connected_clusters
(**kwargs)showdendro
([cores_idx, savefile])Show the clustered dendrogram. Methods Documentation
-
plot_connected_clusters
(**kwargs)[source] [edit on github]¶
-
showdendro
(cores_idx=[], savefile=None)[source] [edit on github]¶ Show the clustered dendrogram. Every color correspond to a different cluster.