fil_finder_2D

class fil_finder.fil_finder_2D(image, header=None, beamwidth=None, ang_scale=None, skel_thresh=None, branch_thresh=None, pad_size=0, skeleton_pad_size=1, flatten_thresh=None, smooth_size=None, size_thresh=None, glob_thresh=None, adapt_thresh=None, distance=None, region_slice=None, mask=None, freq=None, save_name='FilFinder_output')[source]

Bases: BaseInfoMixin

This class acts as an overall wrapper to run the fil-finder algorithm on 2D images and contains visualization and saving capabilities.

Parameters
imagenumpy.ndarray or astropy.io.fits.PrimaryHDU

A 2D array of the data to be analyzed. If a FITS HDU is passed, the header is automatically loaded.

headerFITS header, optional

The header from fits file containing the data. If no header is provided, and it could not be loaded from image, all results will be returned in pixel units.

beamwidthfloat or astropy.units.Quantity, optional

The FWHM beamwidth with an appropriately attached unit. By default, the beam is read from a provided header. If the beam cannot be read from the header, or a header is not provided, this input must be given. If a float is given, it is assumed to be in pixel units.

ang_scaleQuantity, optional

Give the angular to pixel units conversion. If none is given, it will be read from the header. The units must be a valid angular unit.

skel_threshfloat, optional

Given in pixel units.Below this cut off, skeletons with less pixels will be deleted. The default value is 0.3 pc converted to pixels.

branch_threshfloat, optional

Any branches shorter than this length (in pixels) will be labeled as extraneous and pruned off. The default value is 3 times the FWHM beamwidth.

pad_sizeint, optional

The size of the pad (in pixels) used to pad the individual filament arrays. By default this is disabled. If the data continue up to the edge of the image, padding should not be enabled as this causes deviations in the mask around the edges!

skeleton_pad_sizeint, optional

Number of pixels to pad the individual skeleton arrays by. For the skeleton to graph conversion, the pad must always be greater then 0. The amount of padding can slightly effect the extent of the radial intensity profile..

flatten_threshint, optional

The percentile of the data (0-100) to set the normalization of the arctan transform. By default, a log-normal distribution is fit and the threshold is set to \(\mu + 2\sigma\). If the data contains regions of a much higher intensity than the mean, it is recommended this be set >95 percentile.

smooth_sizeint, optional

The patch size (in pixels) used to smooth the flatten image before adaptive thresholding is performed. Smoothing is necessary to ensure the extraneous branches on the skeletons is minimized. If None, the patch size is set to ~0.05 pc. This ensures the large scale structure is not affected while smoothing extraneous pixels off the edges.

size_threshint, optional

This sets the lower threshold on the size of objects found in the adaptive thresholding. If None, the value is set at \(5\pi (0.1 ext(pc))^2\) which is the area of the minimum dimensions expected for a filament. Any region smaller than this threshold may be safely labeled as an artifact of the thresholding.

glob_threshfloat, optional

This is the percentile of the data to mask off. All intensities below are cut off from being included in the filamentary structure.

adapt_threshint, optional

This is the size in pixels of the patch used in the adaptive thresholding. Bright structure is not very sensitive to the choice of patch size, but faint structure is very sensitive. If None, the patch size is set to twice the width of a typical filament (~0.2 pc). As the width of filaments is somewhat ubiquitous, this patch size generally segments all filamentary structure in a given image.

distancefloat, optional

The distance to the region being examined (in pc). If None, the analysis is carried out in pixel and angular units. In this case, the physical priors used in other optional parameters is meaningless and each must be specified initially.

region_slicelist, optional

This gives the option to examine a specific region in the given image. The expected input is [xmin,xmax,ymin,max].

masknumpy.ndarray, optional

A pre-made, boolean mask may be supplied to skip the segmentation process. The algorithm will skeletonize and run the analysis portions only.

freqfloat, optional

Deprecated. Has no effect.

save_namestr, optional

Sets the prefix name that is used for output files. Can be overridden in save_fits and save_table. Default is “FilFinder_output”.

Examples

>>> from fil_finder import fil_finder_2D
>>> from astropy.io import fits
>>> import astropy.units as u
>>> img,hdr = fits.open("twod.fits")[0] 
>>> filfind = fil_finder_2D(img, header=hdr, beamwidth=15*u.arcsec, distance=170*u.pc, save_name='twod_filaments') 
>>> filfind.run(verbose=False) 

Attributes Summary

flat_img_nopad

header

FITS Header.

image_nopad

mask_nopad

medial_axis_distance_nopad

pad_size

skeleton_longpath_nopad

skeleton_nopad

skeleton_pad_size

wcs

WCS Object.

Methods Summary

analyze_skeletons([prune_criteria, ...])

This function wraps most of the skeleton analysis.

compute_filament_brightness()

Returns the median brightness along the skeleton of the filament.

create_mask([glob_thresh, adapt_thresh, ...])

This runs the complete segmentation process and returns a mask of the filaments found.

exec_rht([radius, ntheta, ...])

Implements the Rolling Hough Transform (Clark et al., 2014).

filament_model([max_radius, use_nopad])

Returns a model of the diffuse filamentary network based on the radial profiles.

find_covering_fraction([max_radius])

Compute the fraction of the intensity in the image contained in the filamentary structure.

find_widths([fit_model, try_nonparam, ...])

The final step of the algorithm is to find the widths of each of the skeletons.

medskel([verbose, save_png])

This function performs the medial axis transform (skeletonization) on the mask.

run([verbose, save_name, save_png, table_type])

The whole algorithm in one easy step.

save_fits([save_name, stamps, filename, ...])

This function saves the mask and the skeleton array as FITS files.

save_table([table_type, path, save_name, ...])

The results of the algorithm are saved as an Astropy table in a 'csv', 'fits' or latex format.

Attributes Documentation

flat_img_nopad
header
image_nopad
mask_nopad
medial_axis_distance_nopad
pad_size
skeleton_longpath_nopad
skeleton_nopad
skeleton_pad_size
wcs

Methods Documentation

analyze_skeletons(prune_criteria='all', relintens_thresh=0.2, nbeam_lengths=5, branch_nbeam_lengths=3, skel_thresh=None, branch_thresh=None, verbose=False, save_png=False)[source]

This function wraps most of the skeleton analysis. Several steps are completed here:

  • isolatefilaments is run to separate each skeleton into its own array. If the skeletons are under the threshold set by self.size_thresh, the region is removed. An updated mask is also returned.

  • pix_identify classifies each of the pixels in a skeleton as a body, end, or intersection point. See the documentation on find_filpix for a complete explanation. The function labels the branches and intersections of each skeletons.

  • init_lengths finds the length of each branch in each skeleton and also returns the coordinates of each of these branches for use in the graph representation.

  • pre_graph turns the skeleton structures into a graphing format compatible with networkx. Hubs in the graph are the intersections and end points, labeled as letters and numbers respectively. Edges define the connectivity of the hubs and they are weighted by their length.

  • longest_path utilizes networkx.shortest_path_length to find the overall length of each of the filaments. The returned path is the longest path through the skeleton. If loops exist in the skeleton, the longest path is chosen (this shortest path algorithm fails when used on loops).

  • extremum_pts returns the locations of the longest path’s extent so its performance can be evaluated.

  • final_lengths takes the path returned from longest_path and calculates the overall length of the filament. This step also acts as to prune the skeletons.

  • final_analysis combines the outputs and returns the results for further analysis.

Parameters
verbosebool, optional

Enables plotting.

prune_criteria{‘all’, ‘intensity’, ‘length’}, optional

Choose the property to base pruning on. ‘all’ requires that the branch fails to satisfy the length and relative intensity checks.

relintens_threshfloat, optional

Relative intensity threshold for pruning. Sets the importance a branch must have in intensity relative to all other branches in the skeleton. Must be between (0.0, 1.0].

nbeam_lengthsfloat or int, optional

Sets the minimum skeleton length based on the number of beam sizes specified.

branch_nbeam_lengthsfloat or int, optional

Sets the minimum branch length based on the number of beam sizes specified.

skel_threshfloat, optional

Manually set the minimum skeleton threshold. Overrides all previous settings.

branch_threshfloat, optional

Manually set the minimum branch length threshold. Overrides all previous settings.

save_pngbool, optional

Saves the plot made in verbose mode. Disabled by default.

Attributes
filament_arrayslist of numpy.ndarray

Contains individual arrays of each skeleton

number_of_filamentsint

The number of individual filaments.

array_offsetslist

A list of coordinates for each filament array.This will be used to recombine the final skeletons into one array.

filament_extentslist

This contains the coordinates of the initial and final position of the skeleton’s extent. It may be used to test the performance of the shortest path algorithm.

lengthslist

Contains the overall lengths of the skeletons

labeled_fil_arrayslist of numpy.ndarray

Contains the final labeled versions of the skeletons.

branch_propertiesdict

The significant branches of the skeletons have their length and number of branches in each skeleton stored here. The keys are: filament_branches, branch_lengths

compute_filament_brightness()[source]

Returns the median brightness along the skeleton of the filament.

Attributes
filament_brightnesslist

Average brightness/intensity over the skeleton pixels for each filament.

create_mask(glob_thresh=None, adapt_thresh=None, smooth_size=None, size_thresh=None, verbose=False, test_mode=False, regrid=True, border_masking=True, zero_border=False, fill_hole_size=None, use_existing_mask=False, save_png=False)[source]

This runs the complete segmentation process and returns a mask of the filaments found. The process is broken into six steps:

  • An arctan tranform is taken to flatten extremely bright regions. Adaptive thresholding is very sensitive to local intensity changes and small, bright objects(ie. cores) will leave patch-sized holes in the mask.

  • The flattened image is smoothed over with a median filter. The size of the patch used here is set to be much smaller than the typical filament width. Smoothing is necessary to minimizing extraneous branches when the medial axis transform is taken.

  • A binary opening is performed using an 8-connected structure element. This is very successful at removing small regions around the edge of the data.

  • Objects smaller than a certain threshold (set to be ~1/10 the area of a small filament) are removed to ensure only regions which are sufficiently large enough to be real structure remain.

The parameters for this function are as previously defined. They are included here for fine-tuning purposes only.

Parameters
smooth_sizeint, optional

See definition in fil_finder_2D inputs.

size_threshint, optional

See definition in fil_finder_2D inputs.

glob_threshfloat, optional

See definition in fil_finder_2D inputs.

adapt_threshint, optional

See definition in fil_finder_2D inputs.

verbosebool, optional

Enables plotting. Default is False.

test_modebool, optional

Plot each masking step. Default is False.

regridbool, optional

Enables the regridding of the image to larger sizes when the patch size for the adaptive thresholding is less than 40 pixels. This decreases fragmentation of regions due to pixellization effects. Default is True.

border_maskingbool, optional

Dilates a mask of the regions along the edge of the image to remove regions dominated by noise. Disabling leads to regions characterized at the image boundaries and should only be used if there is not significant noise at the edges. Default is True.

zero_borderbool, optional

Replaces the NaN border with zeros for the adaptive thresholding. This is useful when emission continues to the edge of the image. Default is False.

fill_hole_sizeint or float, optional

Sets the maximum hole size to fill in the skeletons. If <1, maximum is that proportion of the total number of pixels in skeleton. Otherwise, it sets the maximum number of pixels. Defaults to a square area with length of the beamwidth.

use_existing_maskbool, optional

If mask is already specified, enabling this skips recomputing the mask.

save_pngbool, optional

Saves the plot made in verbose mode. Disabled by default.

Attributes
masknumpy.ndarray

The mask of filaments.

exec_rht(radius=10, ntheta=180, background_percentile=25, branches=False, min_branch_length=3, verbose=False, save_png=False)[source]

Implements the Rolling Hough Transform (Clark et al., 2014). The orientation of each filament is denoted by the mean value of the RHT, which from directional statistics can be defined as: \(\langle\theta \rangle = \frac{1}{2} \tan^{-1}\left(\frac{\Sigma_i w_i\sin2\theta_i}{\Sigma_i w_i\cos2\theta_i}\right)\) where \(w_i\) is the normalized value of the RHT at \(\theta_i\). This definition assumes that \(\Sigma_iw_i=1\). \(\theta\) is defined on \(\left[-\pi/2, \pi/2\right)\). “Curvature” is represented by the IQR confidence interval about the mean, \(\langle\theta \rangle \pm \sin^{-1} \left( u_{\alpha} \sqrt{ \frac{1-\alpha}{2R^2} } \right)\) where \(u_{\alpha}\) is the z-score of the two-tail probability, \(\alpha=\Sigma_i\cos{\left[2w_i\left(\theta_i-\langle\theta\rangle\right)\right]}\) is the estimated weighted second trigonometric moment and \(R^2=\left[\left(\Sigma_iw_i\sin{\theta_i}\right)^2 +\left(\Sigma_iw_i\cos{\theta_i}\right)^2\right]\) is the weighted length of the vector.

These equations can be found in Fisher & Lewis (1983).

Parameters
radiusint

Sets the patch size that the RHT uses.

nthetaint, optional

The number of bins to use for the RHT.

backgroundint, optional

RHT distribution often has a constant background. This sets the percentile to subtract off.

branchesbool, optional

If enabled, runs the RHT on individual branches in the skeleton.

min_branch_lengthint, optional

Sets the minimum pixels a branch must have to calculate the RHT

verbosebool, optional

Enables plotting.

save_pngbool, optional

Saves the plot made in verbose mode. Disabled by default.

References

Clark et al. (2014) Fisher & Lewis (1983)

Attributes
rht_curvaturedict

Contains the median and IQR for each filament.

filament_model(max_radius=25, use_nopad=True)[source]

Returns a model of the diffuse filamentary network based on the radial profiles.

Parameters
max_radiusint, optional

Number of pixels to extend profiles to.

use_nopadbool, optional

Returns the unpadded image size when enabled. Enabled by default.

Returns
model_imagenumpy.ndarray

Array of the model

find_covering_fraction(max_radius=25)[source]

Compute the fraction of the intensity in the image contained in the filamentary structure.

Parameters
max_radiusint, optional

Passed to filament_model

Attributes
covering_fractionfloat

Fraction of the total image intensity contained in the filamentary structure (based on the local, individual fits)

find_widths(fit_model=<function gauss_model>, try_nonparam=True, use_longest_paths=False, verbose=False, save_png=False, **kwargs)[source]

The final step of the algorithm is to find the widths of each of the skeletons. We do this by:

  • A Euclidean Distance Transform is performed on each skeleton. The skeletons are also recombined onto a single array. The individual filament arrays are padded to ensure a proper radial profile is created. If the padded arrays fall outside of the original image, they are trimmed.

  • A user-specified model is fit to each of the radial profiles. The default is a Gaussian with a constant background. This returns the width and central intensity of each filament. The reported widths are the deconvolved FWHM of the gaussian width. For faint or crowded filaments, the fit can fail due to lack of data to fit to. In this case, we estimate the width non-parametrically.

Parameters
fit_modelfunction

Function to fit to the radial profile.

try_nonparambool, optional

If True, uses a non-parametric method to find the properties of the radial profile in cases where the model fails.

use_longest_pathsbool, optional

Optionally use the longest path skeletons for the width fitting. Note that this will disregard all branches off of the longest path.

verbosebool, optional

Enables plotting.

save_pngbool, optional

Saves the plot made in verbose mode. Disabled by default.

Attributes
width_fitsdict

Contains the fit parameters and estimations of the errors from each fit.

total_intensitylist

Sum of the intensity in each filament within 1 FWHM of the skeleton.

medskel(verbose=False, save_png=False)[source]

This function performs the medial axis transform (skeletonization) on the mask. This is essentially a wrapper function of skimage.morphology.medial_axis with the ability to delete narrow regions in the mask.

If the distance transform is returned from the transform, it is used as a pruning step. Regions where the width of a region are far too small (set to >0.01 pc) are deleted. This ensures there no unnecessary connections between filaments.

Parameters
verbosebool, optional

Enables plotting.

save_pngbool, optional

Saves the plot made in verbose mode. Disabled by default.

Attributes
skeletonnumpy.ndarray

The array containing all of the skeletons.

medial_axis_distancenumpy.ndarray

The distance transform used to create the skeletons.

run(verbose=False, save_name=None, save_png=False, table_type='fits')[source]

The whole algorithm in one easy step. Individual parameters have not been included in this batch run. If fine-tuning is needed, it is recommended to run each step individually.

Parameters
verbosebool, optional

Enables the verbose option for each of the steps.

save_namestr, optional

The prefix for the saved file. If None, the name from the header is used.

save_pngbool, optional

Saves the plot made in verbose mode. Disabled by default.

table_typestr, optional

Sets the output type of the table. Supported options are “csv”, “fits” and “latex”.

save_fits(save_name=None, stamps=False, filename=None, model_save=True, **kwargs)[source]

This function saves the mask and the skeleton array as FITS files. Included in the header are the setting used to create them.

Parameters
save_namestr, optional

The prefix for the saved file. If None, the save name specified when fil_finder_2D was first called.

stampsbool, optional

Enables saving of individual stamps

filenamestr, optional

File name of the image used. If None, assumes save_name is the file name.

model_savebool, optional

When enabled, calculates the model using filament_model and saves it in a FITS file.

save_table(table_type='csv', path=None, save_name=None, save_branch_props=True, branch_table_type='hdf5', hdf5_path='data')[source]

The results of the algorithm are saved as an Astropy table in a ‘csv’, ‘fits’ or latex format.

Parameters
table_typestr, optional

Sets the output type of the table. Supported options are “csv”, “fits”, “latex” and “hdf5”.

pathstr, optional

The path where the file should be saved.

save_namestr, optional

The prefix for the saved file. If None, the save name specified when fil_finder_2D was first called.

save_branch_propsbool, optional

When enabled, saves the lists of branch lengths and intensities in a separate file(s). Default is enabled.

branch_table_typestr, optional

Any of the accepted table_types will work here. If using HDF5, just one output file is created with each stored within it.

hdf5_pathstr, optional

Path name within the HDF5 file.

Attributes
dataframeastropy.Table

The dataframe is returned for use with the Analysis class.