Photometry

To measure the Petrosian properties of galaxies, we construct a photometric curve of growth using a set of concentric apertures with varying radii to measure the flux. From this curve of growth, we will measure different properties of the galaxy including the total flux, characteristic radii such as the Petrosian radius and half-light radius, and the concentration index.

To achieve an accurate measurement, we must remove the signature of any background flux and minimize the effects of neighboring sources by masking or removing them.

Loading Example Data

The following data is a cutout of a group of bright galaxies in Abell 2744. The original data was acquired by the Hubble Frontier Fields team via the WFC3 instrument in the F105W filter and can be directly downloaded from the Mikulski Archive for Space Telescopes.

We first use astropy’s CCDData to load the example data and visualize it through matplotlib.

[1]:
from astropy.nddata import CCDData

image = CCDData.read('data/abell_2744_dwarf_galaxy_f105w.fits.gz')
[3]:
from matplotlib import pyplot as plt

plt.rcParams['figure.figsize'] = [10, 10]
plt.rcParams['image.origin'] = 'lower'

vmax = image.data.std() # Use the image std as max and min of all plots
vmin = - vmax

plt.imshow(image.data, vmin=0, vmax=image.data.std())
plt.title("Galaxy in Abell 2744")
plt.xlabel("Pixels")
plt.ylabel("Pixels")
plt.show()
_images/photometry_and_petrosian_4_0.png

Estimate data noise at dark area

In this section, we estimate the noise levels in the image. This estimate is used as a threshold for the detection and segmentation steps. To estimate the background noise, we find a dark region in the image (in this case the corner of the image) and make a cutout using Cutout2D. We then compute the statistics for that cutout region.

[4]:
from astropy.nddata import Cutout2D

# Estimate data noise at dark area
# --------------------------------
noise_cutout_pos = (50, 50)
noise_cutout_size = 50
noise_cutout = Cutout2D(image, noise_cutout_pos, noise_cutout_size)

noise_mean = noise_cutout.data.mean()
noise_sigma = noise_cutout.data.std()
noise_3_sigma = noise_sigma * 3.
noise_8_sigma = noise_sigma * 8.

print(noise_mean, noise_3_sigma, noise_8_sigma)

# Plot image and noise distribution
# ---------------------------------
plt.imshow(noise_cutout.data, vmax=noise_mean+noise_3_sigma, vmin=noise_mean-noise_3_sigma)
plt.title("Dark Patch")
plt.xlabel("Pixels")
plt.ylabel("Pixels")
plt.show()

n, bins, patches = plt.hist(noise_cutout.data.flatten(), bins=35, align='left', color='black')
plt.plot(bins[:-1], n, c='r', linewidth=3)
plt.axvline(noise_mean, label="noise_mean", linestyle="--")

plt.xlabel('Flux Bins [{}]'.format(str(image.unit)))
plt.ylabel('Count')
plt.title('Noise Histogram')
plt.legend()

plt.show()
-0.0023249583 0.004241256276145577 0.011310016736388206
_images/photometry_and_petrosian_7_1.png
_images/photometry_and_petrosian_7_2.png

We can also use the standard deviation from the astropy sigma_clipped_stats function to estimate the detection threshold.

[5]:
from astropy.stats import sigma_clipped_stats
image_mean, image_median, image_stddev = sigma_clipped_stats(image.data, sigma=3)

Catalogs

Before we can perform photometry, we need to first identify the sources in the image. To do this, we use photutils and the wrappers in petrofit.

Make Catalog

We first start by defining the detection threshold and we select this value to be the data standard dev or noise_8_sigma computed in the last section.

[6]:
# Define detect threshold
threshold = noise_8_sigma

Next, we define the parameters of a Gaussian kernel that is used to smooth the image before segmentation. The kernel_size defines the dimensions of the smoothing kernel and fwhm which is used to compute the sigma value of the Gaussian distribution. The fwhm is used to compute the sigma value as follows:

sigma = fwhm * gaussian_fwhm_to_sigma
kernel = Gaussian2DKernel(sigma, x_size=kernel_size, y_size=kernel_size)
[7]:
# Define smoothing kernel
kernel_size = 3
fwhm = 3

We also have to define the number of pixels that make up the smallest object. npixels is the number of connected pixels, each greater than the threshold value, that an object must have to be detected. npixels must be a positive integer.

[8]:
npixels = 4**2

To make a catalog of sources and segmentation maps, we use petrofit’s make_catalog function. The function returns a source catalog (photutils.SourceCatalog), a segmentation image, and, if deblend=True, a deblended segmentation image.

The make_catalog function wraps three steps into one function:

  1. Segmentation:

    To identify sources in the image, we first segment the image. The image is smoothed with a gaussian kernel if kernel_size is provided and clipped at the threshold specified (after smoothing if applicable). The image is then segmented using the petrofit.segmentation.make_segments function, which is a wrapper for photutil’s detect_sources functionality.

  2. Deblending:

    To further distinguish the sources, we use photutils.deblend_sources to deblend the sources into individual galaxies. The contrast parameter is the fraction of the total (blended) source flux that a local peak must have (at any one of the multi-thresholds) to be considered as a separate object and nlevels is the number of multi-thresholding levels to use. If deblend is set to False, None is returend for segm_deblend.

  3. Source Catalog:

    Now that we have deblended the sourecs into individual sources, the next step is to create a source catalog (photutils.SourceCatalog) that contains properties like (xcentroid, ycentroid), eccentricity and area. Note that the deblended map is used to make the source catalog but if deblend is set to False, the segmentation map is used instead.

By setting the plot flag to True, we see plots of the segmentation and deblended segmentation maps. The image is plotted along with a color-coded overplot of the segmentation map (each source is a different color). vmax and vmin can be used the same way as in plt.imshow.

After the computation, the following objects are returned:

  • cat : A catalog of sources.

  • segm : Segmentation map.

  • segm_deblend : Deblended segmentation map.

[9]:
from petrofit.segmentation import make_catalog, plot_segments

cat, segm, segm_deblend = make_catalog(
    image.data,
    threshold,
    deblend=True,
    kernel_size=kernel_size,
    fwhm=fwhm,
    npixels=npixels,
    plot=True, vmax=vmax, vmin=vmin
)
WARNING: AstropyDeprecationWarning: "filter_kernel" was deprecated in version 1.2 and will be removed in a future version. Use argument "kernel" instead. [petrofit.segmentation]
WARNING: AstropyDeprecationWarning: "filter_kernel" was deprecated in version 1.2 and will be removed in a future version. Use argument "kernel" instead. [petrofit.segmentation]
_images/photometry_and_petrosian_17_1.png

To demonstrate the useful information in the catalog, we convert the SourceCatalog to an astropy.table.Table and display the first 10 objects.

[10]:
# Display source properties
print("Num of Targets:", len(cat))

# Convert to table
cat_table = cat.to_table()

cat_table[:10]
Num of Targets: 7
[10]:
QTable length=7
labelxcentroidycentroidsky_centroidbbox_xminbbox_xmaxbbox_yminbbox_ymaxareasemimajor_sigmasemiminor_sigmaorientationeccentricitymin_valuemax_valuelocal_backgroundsegment_fluxsegment_fluxerrkron_fluxkron_fluxerr
pix2pixpixdeg
int64float64float64objectint64int64int64int64float64float64float64float64float64float32float32float64float32float64float64float64
1193.275263063171251.8585802859810874None1851990559.03.77393429499739061.3250245077349381.06934147672940520.93633834769960470.0109302080.0237104840.00.9051196nan3.327252817358581nan
261.062221104906726141.06888923864594None596313814432.01.73388788364399531.248006940521428580.894612828134830.69420874741205940.01103714950.0268458330.00.56128824nan1.8057395975075767nan
3192.40371408536166163.40588574852066None18819715916882.02.27698946379877441.944403614117255853.49395936994620.52037779868623170.0114056420.09297860.02.8309908nan2.415208714789045nan
4108.73654703839398183.93651734361652None10611118018736.01.75721274773469971.2602644632531526-85.337033786891280.6968723802238990.010464810.045437150.00.82629424nan1.0688047050210376nan
5178.7230259661561163.506754312433046None14319934932668.09.402746070783629.149685532616527-41.0305601007049750.230440031253410340.0100385223.67713980.0421.50592nan417.0362275083247nan
6138.5631529969507589.27757468116197None114166641121797.07.467332517608415.47049015909432643.8067229008998140.68067060656792080.008826073.65455580.0243.40125nan232.84364151957925nan
799.9772265773608599.12324178530918None83118811171081.06.6565110371964196.34077465664813378.442590371229760.30432805420926270.0103180160.445806030.056.42685nan60.52884884094449nan

Plotting Segmentation Maps

To plot segmentations, you can use the plot_segments function included in PetroFit as follows:

[11]:
from petrofit.segmentation import plot_segments

plot_segments(segm, image=image.data, vmax=vmax, vmin=vmin)
_images/photometry_and_petrosian_21_0.png

As you can see, the segmentation resulted in the identification sources but the sources at the center were classified as a single object because they have interconnecting pixels that are above the threshold.

Next we use the same functiton to plot the deblended segmentation map, notice how the central sources are now deblended into individual sources:

[12]:
plot_segments(segm_deblend, image=image.data, vmax=vmax, vmin=vmin)
_images/photometry_and_petrosian_23_0.png

We can also plot the background pixels that are not a part of a source’s segmentation footprint using the plot_segment_residual function. We significantly lower the vmax and vmin values so the background pixels become more apparent. This plot can be used to see if the threshold used to segment the image was too high. If the threshold is high, we would notice bright pixels that are part of the source in this plot.

[13]:
from petrofit.segmentation import plot_segment_residual
plot_segment_residual(segm, image.data, vmax=vmax/5)
_images/photometry_and_petrosian_25_0.png

Photometry on Single Source

The purpose of this step is to perform aperture photometry to construct a curve of growth that we can use for the Petrosian measurements.

Source Selection

For this example, we will focus on a single source. We have included a helper function order_cat that will produce a list of indices sorted by a key (default is ‘area’). We use the order_cat function to identify the source of interest and perform photometry on its cutout.

[14]:
from petrofit.photometry import order_cat

# Sort and get the largest object in the catalog
sorted_idx_list = order_cat(cat, key='area', reverse=True)
idx = sorted_idx_list[1] # index 0 is largest
source = cat[idx]  # get source from the catalog

Aperture Radii

To construct the curve of growth, we measure the photometry using circular and/or elliptical apertures of varying, concentric radii that are centered on the source. Before we perform the photometry, we need to provide a list of radii that will be used to construct the circular and elliptical apertures. To achieve this, we have provided a helper function that takes in the max radius in pixels (max_pix) and the number radii (n). The function will return a list of radii by dividing the range (max_pix/n, max_pix) into n equally spaced integers.

[15]:
from petrofit.photometry import make_radius_list

r_list = make_radius_list(
    max_pix=50, # Max pixel to go up to
    n=50 # the number of radii to produce
)

print(repr(r_list))
array([ 1.,  2.,  3.,  4.,  5.,  6.,  7.,  8.,  9., 10., 11., 12., 13.,
       14., 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25., 26.,
       27., 28., 29., 30., 31., 32., 33., 34., 35., 36., 37., 38., 39.,
       40., 41., 42., 43., 44., 45., 46., 47., 48., 49., 50.])

Photometry Calculation

The source_photometry function is used to perform photometry on a single source (photutils.segmentation.catalog.SourceCatalog). In this example, we provide the function with the source object, the raw image (note that this is the 2D array and not CCDData), segmentation from the last step (SegmentationImage) and the list of radii we made using make_radius_list. Given these parameters and the options below the function performs photometry as follows:

  • The positions (max value of source in pixels), elongation and position-angle of the source are determined from the SourceCatalog object. The position will be used to center the apertures and the elongation and position angles will be used as parameters of the elliptical apertures. Each of the radii will be assigned an elliptical aperture with these parameters.

  • If cutout_size is defined, the code will use it to make a cutout of that size with the source centered. If the cutout_size is larger than the image or contains pixels outside the image, those pixels outside of the image are replaced by np.nan.

  • Once the cutout_size is determined, cutouts of the error map (if provided) and image are produced. Before the raw image is cutout, sources that are not the source of interest are masked using the segmentation map. The mask_background option gives us the ability to also mask pixels that are considered to be background pixels because they do not belong to any source’s segmentation map. All masked pixels are replaced by np.nan and are not counted in the returned area array.

  • If bkg_sub is set to true, a 2D plane is used to fit pixels that are below a specified sigma from the mean using the petrofit.fitting.fit_plane function. The sigma sigma value is used to determine noise pixels. Once the pixels above this value are masked, a 2D plane is fit to determine the background. The 2D plane model is then converted into an image and subtracted from the cutout of the target source. sigma_type is used to set how this sigma value will be used. The sigma_type options are 'clip' and 'bound':

    • 'clip' (default): Uses astropy.stats.sigma_clipping.sigma_clip to clip at the provided sigma std value. Note that sigma in this case is the number of stds above the mean.

    • 'bound' : After computing the mean of the image, clip at mean - sigma and mean + sigma. Note that sigma in this case is a value and not the number of stds above the mean.

  • The resulting image (after being noise subtracted if bkg_sub is set to true) is passed to the petrofit.photometry.photometry_step which constructs the apertures and performs photometry.

After calculating the photometry at each radius, three arrays are returned:

  • flux_arr: Photometric sum in aperture.

  • area_arr: Exact area of the aperture.

  • error_arr: if error map is provided, error of measurements.

[16]:
from petrofit.photometry import source_photometry

# Photomerty
flux_arr, area_arr, error_arr = source_photometry(

    # Inputs
    source, # Source (`photutils.segmentation.catalog.SourceCatalog`)
    image.data, # Image as 2D array
    segm_deblend, # Deblended segmentation map of image
    r_list, # list of aperture radii

    # Options
    cutout_size=max(r_list)*2, # Cutout out size, set to double the max radius
    bkg_sub=True, # Subtract background
    sigma=3, sigma_type='clip', # Fit a 2D plane to pixels within 3 sigma of the mean
    plot=True, vmax=vmax, vmin=vmin, # Show plot with max and min defined above
)
plt.show()
6
_images/photometry_and_petrosian_32_1.png
_images/photometry_and_petrosian_32_2.png
_images/photometry_and_petrosian_32_3.png

If source_photometry’s plot option is set to True, four plots are displayed: - The top left plot shows the cutout with the 2D plane background subtraction and the surrounding sources masked (replaced by np.nan). - The top right plot shows the curve of growth with increasing radius. The red lines represent the aperture radius. - The bottom two plots show the source profile sliced at the center of the image in the y and x direction respectively.

Save Photometry Arrays

There are many ways of saving the photometry results. For example adding the photometry to the source catalog with the radii as columns. One simple way to save the results for a single source is to save it as a csv file using AstroPy Table.

[17]:
from astropy.table import Table

t = Table(
    data=[r_list, flux_arr, area_arr, error_arr],
    names=['r_list', 'flux_arr', 'area_arr', 'error_arr'],
)

t.write('temp/abell_2744_galaxy_f105w_photometry.csv', overwrite=True)

t[:10]
[17]:
Table length=10
r_listflux_arrarea_arrerror_arr
float64float64float64float64
1.07.6939692.301498nan
2.027.2857999.205992nan
3.051.55991520.713483nan
4.075.61340136.82397nan
5.097.3551757.537453nan
6.0116.3603382.853932nan
7.0132.708958112.773408nan
8.0146.70659147.295879nan
9.0158.733019186.421347nan
10.0169.027815230.149811nan

Petrosian

Construct Petrosian from Photometry

In this section, we use photometric values stored in flux, aperture area, and radii arrays to construct a Petrosian object. The following inputs are needed as inputs:

  • r_list: Array of radii in pixels.

  • area_list: Array of aperture areas.

  • flux_list : Array of photometric flux values.

These values should represent the curve of growth and can be computed by using the PetroFit photometry tools.

The user can also specify the eta and epsilon values.

  • eta (default=0.2) is the Petrosian value that defines the Petrosian radius.

  • epsilon (default=2) is used to determine the radius of total flux.

    • r_total_flux = r_petrosian * epsilon

[18]:
from petrofit.petrosian import Petrosian

p = Petrosian(r_list, area_arr, flux_arr)

Petrosian Radii

PetroFit uses the curve of growth of a galaxy’s flux to compute its Petrosian properties such as Petrosian radius and concentration index.

Petrosian Radius

The Petrosian radius is defined as the radius at which the Petrosian profile reaches the Eta (eta, default=0.2) value.

[19]:
p.r_petrosian # in pixels
[19]:
15.023004600920185

Petrosian Total Flux Radius

The Petrosian flux or total flux radius is the radius that ideally encloses all the flux of the galaxy. The Petrosian total flux radius is estimated by multiplying r_petrosian with epsilon (default=2).

r_total_flux = r_petrosian * epsilon

We can use the r_total_flux_arcsec function, by passing it a WCS object, to compute the total flux radius in arcsec.

[20]:
p.r_total_flux # pixels
[20]:
30.04600920184037
[21]:
p.r_total_flux_arcsec(image.wcs) # arcsec
[21]:
1.802760217889647

Petrosian Half-Light Radius

The half-light radius contains half of the galaxy’s total flux. To compute the half-light radius, we find the total flux (flux at r_total_flux) and divide it by half to find the “half flux” or “half-light”. We then find the pixel closest to the half-light value and define it as the half-light radius. Please note that interpolation is used between the inputted flux radii to find the radius that best matches the half-light flux.

We can use the r_half_light_arcsec function, by passing it a WCS object, to compute the half-light radius in arcsec.

[22]:
p.r_half_light # pixels
[22]:
6.29125825165033
[23]:
p.r_half_light_arcsec(image.wcs) # arcsec
[23]:
0.3774754251173728

Fraction of Flux Radius

We can compute a radius that contains a specific fraction of the total flux using the fraction_flux_to_r function. For example we can compute the radius that contains 60% of the total flux as follows:

[24]:
p.fraction_flux_to_r(fraction=0.6) # pixels
[24]:
7.921584316863373

Concentration Index

The concentration index is the ratio of two aperture radii that contain a fraction (percent) of the total flux. It is computed as follows

concentration_index = 5 * np.log10( r(fraction_2) / r(fraction_1) )

The default is set to fraction_1 = 0.2 and fraction_2 = 0.8. The concentration_index function returns the r_fraction_1, r_fraction_2 and concentration_index.

In these examples, we comput the default C2080 and C5090 concentration indices for the input galaxy:

[25]:
r_20, r_80, c2080 = p.concentration_index()  # defualt c2080

r_20, r_80, c2080 # Radii in pixels
[25]:
(2.880576115223045, 13.332666533306663, 3.3271883082731417)
[26]:
r_50, r_90, c5090 = p.concentration_index(
    fraction_1=0.5,
    fraction_2=0.9
)

r_50, r_90, c5090 # Radii in pixels
[26]:
(6.29125825165033, 18.64372874574915, 2.358976312863468)

Total Petrosian Flux

We can also use Petrosian to compute the total Petrosian flux, which is defined as the flux at r_total_flux. If the r_total_flux is outside the photometric aperture radii, np.nan is returned.

[27]:
p.total_flux
[27]:
242.7890663863256

For Hubble data, we can use the flux_to_abmag function to convert flux values into mags by providing a header.

[28]:
from petrofit.photometry import flux_to_abmag

flux_to_abmag(p.total_flux, header=image.header)
[28]:
20.305745055741784

Petrosian Plots

Profile Plot

The Petrosian plot shows the Petrosian profile, the eta valued to define the Petrosian radius and the Petrosian radius. The blue points are the data points and the orange curve is the interpolated values.

[29]:
# Plot the Petrosian profile
p.plot()
plt.show()
_images/photometry_and_petrosian_58_0.png

We can overplot the half-light radius and total flux radius by setting plot_r=True. We can also overplot a normalized flux curve of growth by setting plot_normalized_flux=True.

[30]:
# Plot with radii and overplot normalized flux curve of growth
p.plot(plot_r=True, plot_normalized_flux=True)
plt.show()
_images/photometry_and_petrosian_60_0.png

Image Overplot

Another way to visualize the radii is to overplot them over an image. To do this we first plot the image as usual and use the Petrosian.imshow function to overplot the r_half_light, r_total_flux, r_20 and r_80. The Petrosian.imshow requires the center of the apertures and plots the radii in pixels. Since elliptical apertures were used, we also provide the elongation and orientation (theta) of the apertures. We get these values from the source object and use utility functions (get_source_position, get_source_elong, get_source_theta) to extract them:

[31]:
from petrofit.segmentation import get_source_position, get_source_elong, get_source_theta

position = get_source_position(source)
elong = get_source_elong(source)
theta = get_source_theta(source)

p.imshow(position=position, elong=elong, theta=theta, lw=1.25)

plt.imshow(image.data, vmax=vmax, vmin=vmin)

plt.legend()
plt.show()
_images/photometry_and_petrosian_62_0.png

Change eta and epsilon

We can change the Petrosian eta and epsilon values after the Petrosian object has been initiated by setting their respective attributes. After setting the attributes, all other calculations and plots will use the new values.

In this example we copy the Petrosian object and change the eta and epsilon values to see how the radii change. Note how r_half_light and r_total_flux changed. Since epsilon = 1, r_total_flux = r_petrosian

[32]:
from copy import copy

p_copy = copy(p)
p_copy.eta = 0.3
p_copy.epsilon = 1

print('eta =', p_copy.eta)
print('epsilon =', p_copy.epsilon)
print('r_half_light (old vs new) = {:0.2f} vs {:0.2f}'.format(p.r_half_light, p_copy.r_half_light))
print('r_total_flux (old vs new) = {:0.2f} vs {:0.2f}'.format(p.r_total_flux, p_copy.r_total_flux))
p_copy.plot(plot_r=True)
eta = 0.3
epsilon = 1
r_half_light (old vs new) = 6.29 vs 4.52
r_total_flux (old vs new) = 30.05 vs 10.59
_images/photometry_and_petrosian_64_1.png

Petrosian Corrections

The most accurate value for the total flux is given when the epsilon value is a function of the Sersic index. As sersic index has been shown to be correlated with the concentration index, we can use the measured concentration index to determine the correct value of epsilon to return the the radii at which to measure the total flux.

Load Correction Grid

As a convenience, a grid of values for epsilon for different concentration indices and half-light radii is provided such that the radius given by r_total_flux = r_petrosian * epsilon will theoretically return 99% of the flux for the source.

To load a grid, simply pass its path to the PetrosianCorrection class as follows:

[33]:
from petrofit.petrosian import PetrosianCorrection

pc = PetrosianCorrection("data/concentration_index_grid_f105w_60mas.yaml")

Estimate Sersic Index

We can estimate the Sersic index n by passing the estimate_n function the uncorrected r_half_light and C2080.

[34]:
estimated_n = pc.estimate_n(
    p.r_half_light,
    p.concentration_index()[-1]
)

estimated_n
[34]:
2.0291863731835336

Estimate Corrected Epsilon

epsilon is the factor we multiply the Petrosian radius with to find the total flux radius and in turn the half-light radius. Here we estimate the corrected value that will result in a better estimate of the total flux radius and half-light radius by passing the estimate_epsilon function the uncorrected r_half_light and C2080.

[35]:
estimated_epsilon = pc.estimate_epsilon(
    p.r_half_light,
    p.concentration_index()[-1]
)

estimated_epsilon
[35]:
2.434355758894113

Corrected Petrosian

Now that we have a better estimate of epsilon, we override the default epsilon value. We recompute the half-light radius and expect it to match the r_eff parameter above (r_eff is the half-light radius). To save the original Petrosian object, we make a new instance with the corrected epsilon:

[36]:
p_corrected = Petrosian(
    p.r_list,
    p.area_list,
    p.flux_list,
    epsilon=estimated_epsilon,
)

Plot the corrected Petrosian radii:

[37]:
p_corrected.plot(plot_r=True, plot_normalized_flux=True)
_images/photometry_and_petrosian_75_0.png

We use the same position, elong, and theta values from the Image Overplot section to overplot corrected and uncorrected (in black) radii over the image.

[38]:
# Uncorrected (Black)
p.imshow(position=position, elong=elong, theta=theta, color='black')

# Corrected
p_corrected.imshow(position=position, elong=elong, theta=theta)

plt.imshow(image.data, vmax=vmax, vmin=vmin)
plt.legend()
plt.show()
_images/photometry_and_petrosian_77_0.png

Corrected Petrosian Total Flux

[39]:
print("Uncorrected Flux = {}".format(p.total_flux * image.unit))
print("Corrected Flux = {}".format(p_corrected.total_flux * image.unit))
Uncorrected Flux = 242.7890663863256 electron / s
Corrected Flux = 245.79963808323313 electron / s