FilFinder2D

class fil_finder.FilFinder2D(image, header=None, beamwidth=None, ang_scale=None, distance=None, mask=None, save_name='FilFinder_output', capture_pre_recombine_masks=False)[source]

Bases: BaseInfoMixin

Extract and analyze filamentary structure from a 2D image.

Parameters
imagendarray or 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.

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.

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.

capture_pre_recombine_masksbool, optional

If True, will save the pre-recombine_skeletons() mask objects and corners and expose them as attributes. Default is False.

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 FilFinder2D
>>> from astropy.io import fits
>>> import astropy.units as u
>>> hdu = fits.open("twod.fits")[0] 
>>> filfind = FilFinder2D(hdu, beamwidth=15*u.arcsec, distance=170*u.pc, save_name='twod_filaments') 
>>> filfind.preprocess_image(verbose=False) 
>>> filfind.create_mask(verbose=False) 
>>> filfind.medskel(verbose=False) 
>>> filfind.analyze_skeletons(verbose=False) 
>>> filfind.exec_rht(verbose=False) 
>>> filfind.find_widths(verbose=False) 
>>> fil_table = filfind.output_table(verbose=False) 
>>> branch_table = filfind.branch_tables(verbose=False) 
>>> filfind.save_fits() 
>>> filfind.save_stamp_fits() 

Attributes Summary

curvature

Returns the orientations of the filament longest paths computed with exec_rht with branches=False.

curvature_branches

Returns the orientations of the filament longest paths computed with exec_rht with branches=False.

end_pts

End pixels for each filament.

flatten_threshold

Threshold value used in the arctan transform.

intersec_pts

Intersection pixels for each filament.

orientation

Returns the orientations of the filament longest paths computed with exec_rht with branches=False.

orientation_branches

Returns the orientations of the filament longest paths computed with exec_rht with branches=False.

pre_recombine_mask_corners

Returns the pre-recombine skeletons mask corners.

pre_recombine_mask_objs

Returns the pre recombine_skeletons() mask objects.

Methods Summary

analyze_skeletons([prune_criteria, ...])

Prune skeleton structure and calculate the branch and longest-path lengths.

branch_lengths([unit])

Return the length of all branches in all filaments.

branch_tables([include_rht])

Return the branch properties of each filament.

covering_fraction([max_radius, ...])

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

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, bkg_subtract, ...])

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

filament_positions([world_coord])

Return the median pixel or world positions of the filaments.

find_widths([max_dist, pad_to_distance, ...])

Create an average radial profile for each filaments and fit a given model.

lengths([unit])

Return longest path lengths of the filaments.

median_brightness()

Returns the median brightness along the skeleton of the filament.

medskel([verbose, save_png])

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

output_table([xunit, world_coord])

Return the analysis results as an astropy table.

preprocess_image([skip_flatten, flatten_percent])

Preprocess and flatten the image before running the masking routine.

ridge_profiles()

Return the image values along the longest path of the skeleton.

save_fits([save_name, ...])

Save the mask and the skeleton array as FITS files.

save_stamp_fits([save_name, pad_size, ...])

Save stamps of each filament image, skeleton, longest-path skeleton, and the model image.

total_intensity([bkg_subtract, bkg_mod_index])

Return the sum of all pixels within the FWHM of the filament.

width_fits([xunit])

Return an Table of the width fit parameters, uncertainties, and whether a flag was raised for a bad fit.

widths([unit])

Fitted FWHM of the filaments and their uncertainties.

Attributes Documentation

curvature

Returns the orientations of the filament longest paths computed with exec_rht with branches=False.

curvature_branches

Returns the orientations of the filament longest paths computed with exec_rht with branches=False.

end_pts

End pixels for each filament.

flatten_threshold

Threshold value used in the arctan transform.

intersec_pts

Intersection pixels for each filament.

orientation

Returns the orientations of the filament longest paths computed with exec_rht with branches=False.

orientation_branches

Returns the orientations of the filament longest paths computed with exec_rht with branches=False.

pre_recombine_mask_corners

Returns the pre-recombine skeletons mask corners. These corners are needed if you want to use FilFinder2D.pre_recombine_mask_objs.

pre_recombine_mask_objs

Returns the pre recombine_skeletons() mask objects. These objects will only be captured if capture_pre_recombine_masks=True and FilFinder2D.create_mask has been run with use_existing_mask=False. Otherwise will return None. This is useful if there are multiple filamentary objects in the image and you want to extract the masks for individual filaments.

Methods Documentation

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

Prune skeleton structure and calculate the branch and longest-path lengths. See skeleton_analysis.

Parameters
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

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.

max_prune_iterint, optional

Maximum number of pruning iterations to apply.

verbosebool, optional

Enables plotting.

save_pngbool, optional

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

save_namestr, optional

Prefix for the saved plots.

branch_lengths(unit=Unit('pix'))[source]

Return the length of all branches in all filaments.

Parameters
unitUnit, optional

Pixel, angular, or physical unit to convert to.

branch_tables(include_rht=False)[source]

Return the branch properties of each filament. If the RHT was run on individual branches (branches=True in exec_rht), the orientation and curvature of each branch can be included in the saved table.

A table will be returned for each filament in order of the filaments in filaments.

Parameters
include_rhtbool, optional

Include RHT orientation and curvature if exec_rht is run with branches=True.

Returns
tableslist

List of Table for each filament.

covering_fraction(max_radius=None, bkg_subtract=True, bkg_mod_index=2)[source]

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

Parameters
max_radiusQuantity, optional

Number of pixels to extend profiles to. If None is given, each filament model is computed to 3 * FWHM.

bkg_subtractbool, optional

Subtract off the fitted background level.

bkg_mod_indexint, optional

Indicate which element in Filament2D.radprof_params is the background level. Defaults to 2 for the Gaussian with background model.

Returns
covering_fractionfloat

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

create_mask(glob_thresh=None, adapt_thresh=None, smooth_size=None, size_thresh=None, verbose=False, test_mode=False, regrid=True, border_masking=True, border_kwargs={'eros_iter': 15, 'filt_width': <Quantity 25. pix>, 'size': <Quantity 50. pix2>}, 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

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.

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.

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=<Quantity 10. pix>, ntheta=180, background_percentile=25, branches=False, min_branch_length=<Quantity 3. pix>, verbose=False, save_png=False, save_name=None)[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.

save_namestr, optional

Prefix for the saved plots.

References

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

Attributes
rht_curvaturedict

Contains the median and IQR for each filament.

filament_model(max_radius=None, bkg_subtract=True, bkg_mod_index=2)[source]

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

Parameters
max_radiusQuantity, optional

Number of pixels to extend profiles to. If None is given, each filament model is computed to 3 * FWHM.

bkg_subtractbool, optional

Subtract off the fitted background level.

bkg_mod_indexint, optional

Indicate which element in Filament2D.radprof_params is the background level. Defaults to 2 for the Gaussian with background model.

Returns
model_imagendarray

Array of the model

filament_positions(world_coord=False)[source]

Return the median pixel or world positions of the filaments.

Parameters
world_coordbool, optional

Return the world coordinates, defined by the WCS information. If no WCS information is given, the output stays in pixel units.

Returns
filament positionslist of tuples

The median positions of each filament.

find_widths(max_dist=<Quantity 10. pix>, pad_to_distance=<Quantity 0. pix>, fit_model='gaussian_bkg', fitter=None, try_nonparam=True, use_longest_path=False, add_width_to_length=True, deconvolve_width=True, fwhm_function=None, chisq_max=10.0, verbose=False, save_png=False, save_name=None, xunit=Unit("pix"), **kwargs)[source]

Create an average radial profile for each filaments and fit a given model. See width_analysis.

  • Radial profiles are created from a Euclidean Distance Transform on the skeleton.

  • A user-specified model is fit to each of the radial profiles. The default model is a Gaussian with a constant background (‘gaussian_bkg’). Other built-in models include a Gaussian with no background (‘gaussian_nobkg’) or a non-parametric estimate (‘nonparam’). Any 1D astropy model (or compound model) can be passed for fitting.

Parameters
imageQuantity or ndarray

The image from which the filament was extracted.

all_skeleton_arraynp.ndarray

An array with the skeletons of other filaments. This is used to avoid double-counting pixels in the radial profiles in nearby filaments.

max_distQuantity, optional

Largest radius around the skeleton to create the profile from. This can be given in physical, angular, or physical units.

pad_to_distanceQuantity, optional

Force all pixels within this distance to be kept, even if a pixel is closer to another skeleton, as given in all_skeleton_array.

fit_modelstr or Fittable1DModel, optional

The model to fit to the profile. Built-in models include ‘gaussian_bkg’ for a Gaussian with a constant background, ‘gaussian_nobkg’ for just a Gaussian, ‘nonparam’ for the non-parametric estimator. Defaults to ‘gaussian_bkg’.

fitterFitter, optional

One of the astropy fitting classes. Defaults to a Levenberg-Marquardt fitter.

try_nonparambool, optional

If the chosen model fit fails, fall back to a non-parametric estimate.

use_longest_pathbool, optional

Only fit profile to the longest path skeleton. Disabled by default.

add_width_to_lengthbool, optional

Add the FWHM to the filament length. This accounts for the expected shortening in the medial axis transform. Enabled by default.

deconvolve_widthbool, optional

Deconvolve the beam width from the FWHM. Enabled by default.

fwhm_functionfunction, optional

Convert the width parameter to the FWHM. Must take the fit model as an argument and return the FWHM and its uncertainty. If no function is given, the Gaussian FWHM is used.

chisq_maxfloat, optional

Enable the fail flag if the reduced chi-squared value is above this limit.

verbosebool, optional

Enables plotting.

save_pngbool, optional

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

save_namestr, optional

Prefix for the saved plots.

xunitUnit, optional

Pixel, angular, or physical unit to convert to in the plot.

kwargsPassed to radial_profile.
lengths(unit=Unit('pix'))[source]

Return longest path lengths of the filaments.

Parameters
unitUnit, optional

Pixel, angular, or physical unit to convert to.

median_brightness()[source]

Returns the median brightness along the skeleton of the filament.

Returns
filament_brightnesslist

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

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.

output_table(xunit=Unit('pix'), world_coord=False, **kwargs)[source]

Return the analysis results as an astropy table.

If exec_rht was run on the whole skeleton, the orientation and curvature will be included in the table. If the RHT was run on individual branches, use save_branch_tables with include_rht=True to save the curvature and orientations.

Parameters
xunitUnit, optional

Unit for spatial properties. Defaults to pixel units.

world_coordbool, optional

Return the median filament position in world coordinates.

kwargsPassed to total_intensity.
preprocess_image(skip_flatten=False, flatten_percent=None)[source]

Preprocess and flatten the image before running the masking routine.

Parameters
skip_flattenbool, optional

Skip the flattening step and use the original image to construct the mask. Default is False.

flatten_percentint, 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.

ridge_profiles()[source]

Return the image values along the longest path of the skeleton. See ridge_profile.

Returns
ridgeslist

List of the ridge values for each filament.

save_fits(save_name=None, save_longpath_skeletons=True, save_model=True, model_kwargs={}, **kwargs)[source]

Save the mask and the skeleton array as FITS files. The header includes the settings used to create them.

The mask, skeleton, (optional) longest skeletons, and (optional) model are included in the outputted file. The skeletons are labeled to match their order in filaments.

Parameters
save_namestr, optional

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

save_longpath_skeletonsbool, optional

Save a FITS extension with the longest path skeleton array. Default is True. Requires analyze_skeletons to be run.

save_modelbool, optional

Save a FITS extension with the longest path skeleton array. Default is True. Requires find_widths to be run.

model_kwargsdict, optional

Passed to filament_model.

kwargsPassed to writeto.
save_stamp_fits(save_name=None, pad_size=<Quantity 20. pix>, model_kwargs={}, **kwargs)[source]

Save stamps of each filament image, skeleton, longest-path skeleton, and the model image.

A suffix of “stamp_{num}” is added to each file, where the number is is the order in the list of filaments.

Parameters
save_namestr, optional

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

stampsbool, optional

Enables saving of individual stamps

model_kwargsdict, optional

Passed to filament_model.

kwargsPassed to writeto.
total_intensity(bkg_subtract=False, bkg_mod_index=2)[source]

Return the sum of all pixels within the FWHM of the filament.

Warning

fil_finder_2D multiplied the total intensity by the angular size of a pixel. This function is just the sum of pixel values. Unit conversions can be applied on the output if needed.

Parameters
bkg_subtractbool, optional

Subtract off the fitted background level.

bkg_mod_indexint, optional

Indicate which element in Filament2D.radprof_params is the background level. Defaults to 2 for the Gaussian with background model.

Returns
total_intensityQuantity

Array of the total intensities for the filament.

width_fits(xunit=Unit('pix'))[source]

Return an Table of the width fit parameters, uncertainties, and whether a flag was raised for a bad fit.

Parameters
xunitUnit, optional

Pixel, angular, or physical unit to convert to.

Returns
tabTable

Table with width fit results.

widths(unit=Unit('pix'))[source]

Fitted FWHM of the filaments and their uncertainties.

Parameters
unitQuantity, optional

The output unit for the FWHM. Default is in pixel units.