title

Photometry of a JWST/NIRCam extragalactic field with Photutils#

Author: Larry Bradley, Space Telescope Science Institute
Data last modified: 2025-02-27

Tutorial Overview#

This tutorial will demonstrate how to use Photutils to perform multiband NIRCam photometry on an extragalactic field.

Photutils is a Python library that provides commonly-used tools and key functionality for detecting and performing photometry of astronomical sources. Tools are provided for

  • background estimation

  • star finding

  • source detection and extraction

  • aperture photometry

  • PSF photometry

  • image segmentation

  • centroids

  • radial profiles

  • elliptical isophote fitting

It is a coordinated package of Astropy and integrates well with other Astropy packages, making it a powerful tool for astronomical image analysis.

Outline:#

  1. Load Datasets

  2. Background Subtraction

  3. Aperture Photometry

  4. Creating a multiband segmentation-based catalog

    1. Creating a Detection Image

    2. Source Extraction via Image Segmentation

    3. Source Deblending

    4. Modifying a Segmentation Image

    5. SourceFinder convenience class

    6. Source Measurements: Source Catalog

    7. Creating a Multiband Catalog

    8. Catalog Exploration

Import packages#

import os
import subprocess
from astroquery.mast import Observations
import astropy.units as u
import matplotlib as mpl
import matplotlib.pyplot as plt
import numpy as np
from astropy.convolution import convolve
from astropy.io import fits
from astropy.table import hstack
from astropy.visualization import SimpleNorm, simple_norm
from astropy.wcs import WCS
from photutils.aperture import (ApertureStats, CircularAperture,
                                SkyCircularAperture, aperture_photometry, 
                                aperture_to_region)
from photutils.background import (Background2D, BiweightLocationBackground,
                                  BiweightScaleBackgroundRMS)
from photutils.segmentation import (SourceCatalog,
                                    deblend_sources, detect_sources,
                                    make_2dgaussian_kernel)
from photutils.utils import CutoutImage
from regions import Regions
from scipy.ndimage import grey_dilation

# Change some default plotting parameters
mpl.rcParams['figure.constrained_layout.use'] = True
mpl.rcParams['figure.max_open_warning'] = 40
mpl.rcParams['image.interpolation'] = 'nearest'
mpl.rcParams['image.origin'] = 'lower'

# Enable interactive plots
# %matplotlib ipympl
def download_files(files_to_download):
    """Download a list of files from MAST.

    Parameters
    ----------
    files_to_download : list
        List of filenames
    """
    for file in files_to_download:
        # Check if the file already exists in the current working directory
        if os.path.exists(file):
            print(f"File {file} already exists. Skipping download.")
            continue
        cal_uri = f'mast:HLSP/ceers/nircam/{file}.gz'
        Observations.download_file(cal_uri)
        command = ['gunzip', file, '.gz']
        subprocess.run(command, capture_output=True, text=True)

Load Datasets#

Data Overview#

We will use 3 JWST/NIRCam images in the F115W, F200W, and F356W filters from the Cosmic Evolution Early Release Science Survey (CEERS) survey.

These are fully-calibration High-level Science Products (HLSPs) produced by the CEERS team of their NIRCam 2 pointing within the Extended Groth Strip (EGS) HST legacy field. The full CEERS data includes 7 NIRCam filters (F115W, F150W, F200W, F277W, F356W, F410M, F444W) plus MIRI and NIRSpec data. The NIRCam data were taken in parallel to prime NIRSpec observations.

The images are pixel aligned and have a pixel scale of 30 mas/pixel.

fn1 = 'hlsp_ceers_jwst_nircam_nircam2_f115w_v0.5_i2d.fits'
fn2 = 'hlsp_ceers_jwst_nircam_nircam2_f200w_v0.5_i2d.fits'
fn3 = 'hlsp_ceers_jwst_nircam_nircam2_f356w_v0.5_i2d.fits'
download_files([fn1, fn2, fn3])
Downloading URL https://mast.stsci.edu/api/v0.1/Download/file?uri=mast:HLSP/ceers/nircam/hlsp_ceers_jwst_nircam_nircam2_f115w_v0.5_i2d.fits.gz to hlsp_ceers_jwst_nircam_nircam2_f115w_v0.5_i2d.fits.gz ...
 [Done]
Downloading URL https://mast.stsci.edu/api/v0.1/Download/file?uri=mast:HLSP/ceers/nircam/hlsp_ceers_jwst_nircam_nircam2_f200w_v0.5_i2d.fits.gz to hlsp_ceers_jwst_nircam_nircam2_f200w_v0.5_i2d.fits.gz ...
 [Done]
Downloading URL https://mast.stsci.edu/api/v0.1/Download/file?uri=mast:HLSP/ceers/nircam/hlsp_ceers_jwst_nircam_nircam2_f356w_v0.5_i2d.fits.gz to hlsp_ceers_jwst_nircam_nircam2_f356w_v0.5_i2d.fits.gz ...
 [Done]

Let’s see what’s in these files using the fits.info() function.

fits.info(fn1)
Filename: hlsp_ceers_jwst_nircam_nircam2_f115w_v0.5_i2d.fits
No.    Name      Ver    Type      Cards   Dimensions   Format
  0  PRIMARY       1 PrimaryHDU     374   ()      
  1  SCI_BKSUB     1 ImageHDU        76   (11000, 6450)   float64   
  2  SCI           1 ImageHDU        75   (11000, 6450)   float32   
  3  ERR           1 ImageHDU        10   (11000, 6450)   float32   
  4  CON           1 ImageHDU         9   (11000, 6450)   int32   
  5  WHT           1 ImageHDU         9   (11000, 6450)   float32   
  6  VAR_POISSON    1 ImageHDU         9   (11000, 6450)   float32   
  7  VAR_RNOISE    1 ImageHDU         9   (11000, 6450)   float32   
  8  VAR_FLAT      1 ImageHDU         9   (11000, 6450)   float32   
  9  BKGD          1 ImageHDU        39   (11000, 6450)   float64   
 10  BKGMASK       1 ImageHDU        39   (11000, 6450)   float32   
 11  HDRTAB        1 BinTableHDU    816   96R x 403C   [23A, 5A, 3A, 48A, 6A, 13A, 5A, 5A, 7A, 10A, 4A, L, D, D, D, D, 4A, 18A, 57A, 22A, 3A, D, 20A, D, 10A, 12A, 23A, 23A, 26A, 11A, 5A, 3A, 3A, 2A, 1A, 2A, 1A, L, 12A, 23A, 2A, 26A, 20A, 27A, 10A, K, L, L, L, L, 5A, 1A, 5A, D, D, D, D, D, D, D, D, D, D, 6A, 5A, 1A, 5A, 5A, 5A, L, D, D, D, D, D, D, D, D, D, D, D, D, 4A, D, D, D, D, D, D, D, D, D, K, 20A, 9A, D, D, D, D, D, D, D, D, D, 7A, D, D, K, K, D, D, K, K, D, D, K, K, K, K, K, D, D, D, D, D, D, D, D, K, K, L, L, K, K, D, D, D, D, D, D, D, 4A, K, K, K, K, K, K, D, D, D, D, 4A, D, D, K, D, K, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, 7A, 9A, D, D, D, D, D, D, D, D, D, D, D, D, D, 10A, 11A, D, D, D, D, D, D, D, D, D, D, D, D, K, K, D, 4A, K, K, K, D, 4A, K, K, K, D, 4A, K, D, D, K, 27A, 27A, 10A, D, D, D, D, D, D, D, 9A, 27A, D, D, D, D, D, D, D, 8A, 14A, 33A, D, D, 3A, 3A, D, 33A, 3A, 39A, D, D, 41A, 33A, 3A, 3A, 3A, 3A, 3A, D, 33A, 3A, 3A, 3A, D, D, 38A, 33A, 3A, 3A, D, 35A, 35A, D, 38A, D, 3A, D, D, D, D, 39A, D, D, D, 3A, D, 38A, D, 40A, 37A, D, D, D, 3A, D, D, D, D, D, 8A, D, D, D, D, D, 8A, 8A, D, D, D, D, D, 8A, D, 7A, 7A, D, D, 7A, 8A, D, D, 8A, D, D, D, 7A, D, 8A, 8A, 8A, 8A, 7A, D, D, 8A, 7A, D, D, D, D, 8A, D, 7A, D, D, D, 5A, D, L, 6A, D, D, D, D, 4A, D, D, D, K, D, D, D, D, D, D, 12A, 12A, D, 3A, 3A, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, 121A, D, D, D, D, D, D, K, D, D, D, D]   
 12  ASDF          1 BinTableHDU     11   1R x 1C   [36489B]   

Note that these files have extra extensions that are not in the pipeline-produced i2d files (e.g., SCI_BKSUB, BKGD, BKGMASK)

Let’s load the NIRCam F115W image and display it.

data = fits.getdata(fn1, ext=2)
snorm = SimpleNorm('sqrt', percent=98)
norm = snorm(data[data != 0])
fig, ax = plt.subplots(figsize=(10, 5))
_ = ax.imshow(data, norm=norm)
../../../_images/ea4cb30fcbe4e33bf355838742f14336971e073d163135ad9239eadc9a35896b.png

There are 2 NIRCam modules (A & B), each with 4 NIRCam short-wavelength (SW) detectors. All 8 SW detectors are 2k x 2k. There are also 2 2k x 2k long wavelength detectors, one covering each module. We see the SW detector gaps in each NIRCam module because these were parallel observations to NIRSpec prime observations. The prime observations had only very small dithers.

Data Units#

The JWST calibration pipeline produces images in units of MJy/sr. The unit is stored in the FITS header in the BUNIT key.

header = fits.getheader(fn1, ext=2)
header['BUNIT']
'MJy/sr'
data_unit = u.Unit(header['BUNIT'])
data_unit
\[\mathrm{\frac{MJy}{sr}}\]
# another way to create the unit
u.MJy / u.sr
\[\mathrm{\frac{MJy}{sr}}\]

The pixel area in steradians is stored in the FITS header in the PIXAR_SR key.

header['PIXAR_SR']
2.11539874851881e-14

Let’s use Astropy units to calculate a scale factor to convert our original data in units of MJy/sr to μJy.

unit = u.uJy
unit_conv = ((1 * data_unit) * (header['PIXAR_SR'] * u.sr)).to(unit)
unit_conv
\[0.021153987 \; \mathrm{\mu Jy}\]

Load and Display the Data#

Let’s load the data. We’ll use only one of the NIRCam modules. We’ll also convert the units from MJy/sr to μJy.

# define a tuple of slice objects to extract only 1 NIRCam module
slc = (slice(1650, 6450), slice(0, 4700))


def load_data(filename):
    data = {}
    with fits.open(filename) as hdul:
        data['filter'] = filename.split('_')[6]
        data['sci_bksub'] = hdul[1].data[slc]
        data['sci'] = hdul[2].data[slc]
        data['err'] = hdul[3].data[slc]
        data['var_rnoise'] = hdul[7].data[slc]
        data['wcs'] = WCS(hdul[1].header)[slc]

        # convert the units in-place
        data['sci_bksub'] *= unit_conv.value
        data['sci_bksub'] <<= unit_conv.unit
        data['sci'] *= unit_conv.value
        data['sci'] <<= unit_conv.unit
        data['err'] *= unit_conv.value
        data['err'] <<= unit_conv.unit
        data['var_rnoise'] *= unit_conv.value ** 2
        data['var_rnoise'] <<= unit_conv.unit ** 2
    return data


f115w = load_data(fn1)
f200w = load_data(fn2)
f356w = load_data(fn3)
data = (f115w, f200w, f356w)
WARNING: FITSFixedWarning: 'datfix' made the change 'Set DATE-BEG to '2022-06-21T22:57:30.824' from MJD-BEG.
Set DATE-AVG to '2022-06-22T00:06:18.287' from MJD-AVG.
Set DATE-END to '2022-06-22T01:14:35.286' from MJD-END'. [astropy.wcs.wcs]
WARNING: FITSFixedWarning: 'obsfix' made the change 'Set OBSGEO-L to   -74.783957 from OBSGEO-[XYZ].
Set OBSGEO-B to   -36.777916 from OBSGEO-[XYZ].
Set OBSGEO-H to 1725421070.363 from OBSGEO-[XYZ]'. [astropy.wcs.wcs]
WARNING: FITSFixedWarning: 'datfix' made the change 'Set DATE-BEG to '2022-06-22T02:33:19.452' from MJD-BEG.
Set DATE-AVG to '2022-06-22T02:59:46.710' from MJD-AVG.
Set DATE-END to '2022-06-22T03:26:17.551' from MJD-END'. [astropy.wcs.wcs]
WARNING: FITSFixedWarning: 'obsfix' made the change 'Set OBSGEO-L to   -74.741833 from OBSGEO-[XYZ].
Set OBSGEO-B to   -36.807444 from OBSGEO-[XYZ].
Set OBSGEO-H to 1725573761.369 from OBSGEO-[XYZ]'. [astropy.wcs.wcs]
WARNING: FITSFixedWarning: 'datfix' made the change 'Set DATE-BEG to '2022-06-22T00:09:16.356' from MJD-BEG.
Set DATE-AVG to '2022-06-22T00:41:54.035' from MJD-AVG.
Set DATE-END to '2022-06-22T01:14:35.286' from MJD-END'. [astropy.wcs.wcs]
WARNING: FITSFixedWarning: 'obsfix' made the change 'Set OBSGEO-L to   -74.770114 from OBSGEO-[XYZ].
Set OBSGEO-B to   -36.787605 from OBSGEO-[XYZ].
Set OBSGEO-H to 1725471424.000 from OBSGEO-[XYZ]'. [astropy.wcs.wcs]

Let’s display the images.

snorm = SimpleNorm('sqrt', percent=98)
norm = snorm(data[2]['sci'].value)

nimg = len(data)
fig, ax = plt.subplots(ncols=nimg, figsize=(12, 4), constrained_layout=True)
for i in range(nimg):
    ax[i].imshow(data[i]['sci'].value, norm=norm)
    ax[i].set_title(data[i]['filter'].upper())
../../../_images/215c8d00e44e36f83d6ce102fcbcfcde9c1d4cd260928ba94b654075f0cd6e29.png

Background Subtraction#

Let’s plot the histograms of the pixel values near zero to see if there is any background pedestal.

The CEERS team already subtracted a constant background value from each image.

fig, ax = plt.subplots(ncols=3, figsize=(10, 4))

for i in range(nimg):
    vals = data[i]['sci'].value
    vals = vals[vals != 0]
    ax[i].hist(vals, bins=150, range=(-0.001, 0.001))
    ax[i].set_title(data[i]['filter'].upper())
    ax[i].axvline(0, color='red', ls='dotted')
../../../_images/df31feb1d0142daf4a603bb8d39f55b370e89100355fc0a8ea9beacbbb7e0a10.png

2D Background Estimation#

To remove the varying residual background, we can use the photutils Background2D class.

This class can be used to estimate a 2D background and background RMS noise in an image.

The background is estimated using (sigma-clipped) statistics in each box of a grid that covers the input data to create a low-resolution background map.

The final background map is calculated by interpolating the low-resolution background map.

Source Mask#

To get the best possible background estimate, we need to first create a mask of the sources in the image.

Let’s create a very rough source mask using a single threshold value. Here, we dilate the mask using a 45x45 pixel square footprint to include source wings in the mask.

mask = f115w['sci'].value > 0.001
mask = grey_dilation(mask, size=45)

fig, ax = plt.subplots()
_ = ax.imshow(mask)
../../../_images/f2f9bf64df402db1cf147806904472db0bc4a6d2b388cf2e7d753f4058c7919b.png

Background2D#

bkg_estimator = BiweightLocationBackground()
bkgrms_estimator = BiweightScaleBackgroundRMS()
coverage_mask = f115w['sci'] == 0
bkg = Background2D(f115w['sci'], box_size=20, filter_size=5, exclude_percentile=90,
                   bkg_estimator=bkg_estimator, bkgrms_estimator=bkgrms_estimator,
                   coverage_mask=coverage_mask, mask=mask)
bkgim = bkg.background
fig, ax = plt.subplots(nrows=2, ncols=2, figsize=(8, 8))
ax = ax.ravel()
ax[0].imshow(f115w['sci'].value, norm=norm)
ax[0].set_title('F115W Original')
ax[1].imshow(bkgim.value, norm=norm)
ax[1].set_title('F115W Background')
ax[2].imshow((f115w['sci'] - bkgim).value, norm=norm)
ax[2].set_title('F115W Background Subtracted')
ax[3].imshow(f115w['sci_bksub'].value, norm=norm)
_ = ax[3].set_title('F115W CEERS Background Subtracted')
../../../_images/58ea1ed967bbe3edcdb2b1e25e00414b2bd18cdeb1ae334984b14b9ccc952be9.png
fits.info(fn1)
Filename: hlsp_ceers_jwst_nircam_nircam2_f115w_v0.5_i2d.fits
No.    Name      Ver    Type      Cards   Dimensions   Format
  0  PRIMARY       1 PrimaryHDU     374   ()      
  1  SCI_BKSUB     1 ImageHDU        76   (11000, 6450)   float64   
  2  SCI           1 ImageHDU        75   (11000, 6450)   float32   
  3  ERR           1 ImageHDU        10   (11000, 6450)   float32   
  4  CON           1 ImageHDU         9   (11000, 6450)   int32   
  5  WHT           1 ImageHDU         9   (11000, 6450)   float32   
  6  VAR_POISSON    1 ImageHDU         9   (11000, 6450)   float32   
  7  VAR_RNOISE    1 ImageHDU         9   (11000, 6450)   float32   
  8  VAR_FLAT      1 ImageHDU         9   (11000, 6450)   float32   
  9  BKGD          1 ImageHDU        39   (11000, 6450)   float64   
 10  BKGMASK       1 ImageHDU        39   (11000, 6450)   float32   
 11  HDRTAB        1 BinTableHDU    816   96R x 403C   [23A, 5A, 3A, 48A, 6A, 13A, 5A, 5A, 7A, 10A, 4A, L, D, D, D, D, 4A, 18A, 57A, 22A, 3A, D, 20A, D, 10A, 12A, 23A, 23A, 26A, 11A, 5A, 3A, 3A, 2A, 1A, 2A, 1A, L, 12A, 23A, 2A, 26A, 20A, 27A, 10A, K, L, L, L, L, 5A, 1A, 5A, D, D, D, D, D, D, D, D, D, D, 6A, 5A, 1A, 5A, 5A, 5A, L, D, D, D, D, D, D, D, D, D, D, D, D, 4A, D, D, D, D, D, D, D, D, D, K, 20A, 9A, D, D, D, D, D, D, D, D, D, 7A, D, D, K, K, D, D, K, K, D, D, K, K, K, K, K, D, D, D, D, D, D, D, D, K, K, L, L, K, K, D, D, D, D, D, D, D, 4A, K, K, K, K, K, K, D, D, D, D, 4A, D, D, K, D, K, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, 7A, 9A, D, D, D, D, D, D, D, D, D, D, D, D, D, 10A, 11A, D, D, D, D, D, D, D, D, D, D, D, D, K, K, D, 4A, K, K, K, D, 4A, K, K, K, D, 4A, K, D, D, K, 27A, 27A, 10A, D, D, D, D, D, D, D, 9A, 27A, D, D, D, D, D, D, D, 8A, 14A, 33A, D, D, 3A, 3A, D, 33A, 3A, 39A, D, D, 41A, 33A, 3A, 3A, 3A, 3A, 3A, D, 33A, 3A, 3A, 3A, D, D, 38A, 33A, 3A, 3A, D, 35A, 35A, D, 38A, D, 3A, D, D, D, D, 39A, D, D, D, 3A, D, 38A, D, 40A, 37A, D, D, D, 3A, D, D, D, D, D, 8A, D, D, D, D, D, 8A, 8A, D, D, D, D, D, 8A, D, 7A, 7A, D, D, 7A, 8A, D, D, 8A, D, D, D, 7A, D, 8A, 8A, 8A, 8A, 7A, D, D, 8A, 7A, D, D, D, D, 8A, D, 7A, D, D, D, 5A, D, L, 6A, D, D, D, D, 4A, D, D, D, K, D, D, D, D, D, D, 12A, 12A, D, 3A, 3A, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, D, 121A, D, D, D, D, D, D, K, D, D, D, D]   
 12  ASDF          1 BinTableHDU     11   1R x 1C   [36489B]   

Let’s load the CEERS source mask and show a comparison.

ceers_bkgmask = fits.getdata(fn1, ext=10)[slc]
fig, ax = plt.subplots(ncols=2, figsize=(10, 5))
ax[0].imshow(mask)
ax[0].set_title('Simple Source Mask')
ax[1].imshow(ceers_bkgmask.astype(bool) ^ coverage_mask)
_ = ax[1].set_title('CEERS Source Mask')
../../../_images/f31c857a992ff73c070c5e3177e6b9fc8ce1fae481b8227d0aa239fb80303dc9.png

Let’s compute the background again, but this time we’ll use the CEERS source mask.

bkg = Background2D(f115w['sci'], box_size=10, filter_size=5,
                   exclude_percentile=90, coverage_mask=coverage_mask, mask=ceers_bkgmask,
                   bkg_estimator=bkg_estimator, bkgrms_estimator=bkgrms_estimator)
bkgim = bkg.background

fig, ax = plt.subplots(nrows=2, ncols=2, figsize=(8, 8))
ax = ax.ravel()
ax[0].imshow(f115w['sci'].value, norm=norm)
ax[0].set_title('F115W Original')
ax[1].imshow(bkgim.value, norm=norm)
ax[1].set_title('F115W Background')
ax[2].imshow((f115w['sci'] - bkgim).value, norm=norm)
ax[2].set_title('F115W Background Subtracted')
ax[3].imshow(f115w['sci_bksub'].value, norm=norm)
_ = ax[3].set_title('F115W CEERS Background Subtracted')
../../../_images/20ab9da0d36f4d07c1867eea9972ac844c70bd6cdfbbdcda0dc67624b1a0cf88.png

See the CEERS NIRCam GitHub repository for the scripts used for making the source mask and performing background subtraction.

We’ll use the CEERS background-subtracted (SCI_BKSUB) images for the rest of this tutorial.

Aperture Photometry#

Photutils provides circular, elliptical, and rectangular aperture shapes (plus annulus versions of each).

Further, there are two types of aperture classes, defined either with pixel or sky (celestial) coordinates.

These are the names of the aperture classes that are defined in pixel coordinates:

  • CircularAperture

  • CircularAnnulus

  • EllipticalAperture

  • EllipticalAnnulus

  • RectangularAperture

  • RectangularAnnulus

The aperture classes defined in celestial coordinates have Sky prepended to their names:

  • SkyCircularAperture

  • SkyCircularAnnulus

  • SkyEllipticalAperture

  • SkyEllipticalAnnulus

  • SkyRectangularAperture

  • SkyRectangularAnnulus

Define an Aperture Object#

First, let’s define a circular aperture at the two given (x, y) pixel positions and radius (in pixels).

Photutils has many tools for detecting sources (e.g., see photutils.detection and photutils.segmentation) and measuring source positions (e.g., photutils.centroids). For this example, we’ll use estimated values based on cursor positions.

xypos = [(1262.9, 2781.4), (1285.9, 2939.3)]
aper = CircularAperture(xypos, r=10)
aper
<CircularAperture([[1262.9, 2781.4],
                 [1285.9, 2939.3]], r=10.0)>
fig, ax = plt.subplots(figsize=(6, 6))
ax.imshow(f115w['sci_bksub'].value, norm=norm)
aper.plot(ax=ax, color='red')
ax.set_xlim(1100, 1400)
_ = ax.set_ylim(2700, 3000)
../../../_images/8cecb6bcd53dc9c96929676a5702e223e67841f3e2bf3fc34dc658b51e3986d0.png

aperture_photometry function#

phot = aperture_photometry(f200w['sci_bksub'], aper, error=f200w['err'], wcs=f200w['wcs'])
phot
QTable length=2
idxcenterycentersky_centeraperture_sumaperture_sum_err
deg,deguJyuJy
int64float64float64SkyCoordfloat64float64
11262.92781.4214.92588919088467,52.9418574199865550.073010656363694750.004254292010347718
21285.92939.3214.92401980802305,52.942563875123470.065536498836604990.0042453372733323195
phot['aperture_mag'] = phot['aperture_sum'].to(u.ABmag)
phot
QTable length=2
idxcenterycentersky_centeraperture_sumaperture_sum_erraperture_mag
deg,deguJyuJymag(AB)
int64float64float64SkyCoordfloat64float64float64
11262.92781.4214.92588919088467,52.9418574199865550.073010656363694750.00425429201034771826.741534368116987
21285.92939.3214.92401980802305,52.942563875123470.065536498836604990.004245337273332319526.858791909236327

Multiple positions and sizes#

Multiple aperture objects can be input, BUT they must all have identical positions.

aper2 = CircularAperture(xypos, r=12)
aperture_photometry(f200w['sci_bksub'], (aper, aper2), error=f200w['err'], wcs=f200w['wcs'])
QTable length=2
idxcenterycentersky_centeraperture_sum_0aperture_sum_err_0aperture_sum_1aperture_sum_err_1
deg,deguJyuJyuJyuJy
int64float64float64SkyCoordfloat64float64float64float64
11262.92781.4214.92588919088467,52.9418574199865550.073010656363694750.0042542920103477180.076873325748865710.005082861205748076
21285.92939.3214.92401980802305,52.942563875123470.065536498836604990.00424533727333231950.065559975208257790.005085923868101762

ApertureStats class#

The ApertureStats calculates various statistics (including the sum) and properties using the pixels within the aperture.

Only one aperture object can be input.

apstats = ApertureStats(f200w['sci_bksub'], aper, error=f200w['err'], wcs=f200w['wcs'])
apstats
<photutils.aperture.stats.ApertureStats>
Length: 2
stats_tbl = apstats.to_table()
stats_tbl
QTable length=2
idxcentroidycentroidsky_centroidsumsum_errsum_aper_areacenter_aper_areaminmaxmeanmedianmodestdmad_stdvarbiweight_locationbiweight_midvariancefwhmsemimajor_sigmasemiminor_sigmaorientationeccentricity
deg,deguJyuJypix2pix2uJyuJyuJyuJyuJyuJyuJyuJy2uJyuJy2pixpixpixdeg
int64float64float64SkyCoordfloat64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64
11263.1826631347332781.8937689404197214.92588145725446,52.941858291419340.073010656363694750.004254292010347718314.1592653589793314.0-0.000288897852823565550.00170842949945191210.00023225154263428050.0001366333865631305-5.4602925579169464e-050.00033367933703977730.000194523166611500331.1134189996730531e-070.00013225202242378064.539571179364758e-087.3339057021657283.3578048236655772.850334712856219-53.7983034798555050.5286040913887978
21285.8849558958282939.0652696702764214.92402241551446,52.942562703398130.065536498836604990.0042453372733323195314.1592653589793313.0-0.0002791923225016960.00220594768406640830.000209272155544103389.157144209544037e-05-0.000143829984801885650.00038731159370063170.000158165049361719561.500102706149232e-078.897633733574341e-053.047323059521578e-086.8239775326976413.20479639965636662.5543405760167427.1656355635916460.603930778084528

The aperture statistics can also be accessed as attributes, e.g.,:

apstats.mean
\[[0.00023225154,~0.00020927216] \; \mathrm{\mu Jy}\]

Note that the returned measurements are for the input aperture position. ApertureStats returns the centroid value of the pixels within the input aperture, but the aperture is not recentered at the measured centroid position. However, you can create a new Aperture object using the measured centroid and then re-run ApertureStats.

Sky-based apertures#

Sky-based apertures require an Astropy SkyCoord position and a size as an Astropy Quantity in angular units.

Sky apertures are not defined completely in sky coordinates. They simply use sky coordinates to define the central position, and the remaining parameters are converted to pixels using the pixel scale of the image at the central position. Projection distortions are not taken into account. They are not defined as apertures on the celestial sphere, but rather are meant to represent aperture shapes on an image. If the apertures were defined completely in sky coordinates, their shapes would not be preserved when converting to or from pixel coordinates.

xypos = [(1262.9, 2781.4), (1285.9, 2939.3)]
skycoord = f200w['wcs'].pixel_to_world(*np.transpose(xypos))
skycoord
<SkyCoord (ICRS): (ra, dec) in deg
    [(214.92588919, 52.94185742), (214.92401981, 52.94256388)]>
sky_aper = SkyCircularAperture(skycoord, r=0.3 * u.arcsec)
sky_aper
<SkyCircularAperture(<SkyCoord (ICRS): (ra, dec) in deg
    [(214.92588919, 52.94185742), (214.92401981, 52.94256388)]>, r=0.3 arcsec)>

When performing aperture photometry with a Sky-based aperture, a WCS transformation must be input.

phot = aperture_photometry(f200w['sci_bksub'], sky_aper, error=f200w['err'], wcs=f200w['wcs'])
phot
QTable length=2
idxcenterycentersky_centeraperture_sumaperture_sum_err
deg,deguJyuJy
int64float64float64SkyCoordfloat64float64
11262.89999999927022781.4000000001824214.92588919088467,52.9418574199865550.07301071229352970.0042543119728190954
21285.8999999972822939.3000000007023214.92401980802305,52.942563875123470.065536564280458050.00424535734794874

Creating a multiband segmentation-based catalog#

Creating a Detection Image#

We are ready to detect sources in our multiband data, but which image do we use? Typically, that depends on your science use case. For example, perhaps you are interested in Lyman-break dropout galaxies, in which case a redder filter or filters is more appropriate.

One other consideration is the varying PSF across our multiband data where the PSF is broader at redder-wavelengths (R ~ λ / D). One method to address this is to PSF match all the multiband data to a common PSF (usually the reddest band with the broadest PSF). A common technique is to generate a set of matching kernels that are then convolved with the images, blurring them all to the same resolution as the reddest band. See Photutils PSF Matching for tools to create matching kernels. We do not PSF match our multiband data in the examples that follow.

In general, one should use a detection image with a high signal-to-noise ratio. A common case is to create a deep detection image by combining the images in several filters. We can do this by creating an inverse-variance weighted combination of images. For this example, we’ll use all three filter images. The weights will be the inverse variance of the read noise (the read noise variance is in the VAR_RNOISE extension).

scis = [filt['sci_bksub'] for filt in data]
whts = [1.0 / filt['var_rnoise'] for filt in data]
detimg = np.multiply(scis, whts).sum(axis=0) / np.sum(whts, axis=0) << unit
scis = None  # free memory
whts = None
/tmp/ipykernel_2196/3398041919.py:3: RuntimeWarning: invalid value encountered in divide
  detimg = np.multiply(scis, whts).sum(axis=0) / np.sum(whts, axis=0) << unit
fig, ax = plt.subplots(figsize=(6, 6))
_ = ax.imshow(detimg.value, norm=norm)
../../../_images/e3b60c55cae78fb04effdc672b097c390051c9c5d1c1ca33f3b866b84d057f74.png

Image Segmentation (Source Extraction)#

Image segmentation is a process of assigning a label to every pixel in an image, such that pixels with the same label are part of the same source. Detected sources must have a minimum number of connected pixels that are each greater than a specified threshold value in an image.

The threshold level is usually defined as a multiple of the background noise (sigma level) above the background. This can be specified either as a per-pixel threshold image or a single threshold value for the whole image.

The image is often filtered before thresholding to smooth the noise and maximize the detectability of objects with a shape similar to the filter kernel.

Let’s convolve the background-subtracted detection image with a 2D Gaussian kernel with a FWHM of 3 pixels.

kernel_fwhm = 3.0
kernel_size = int(kernel_fwhm * 2 + 1)
kernel = make_2dgaussian_kernel(kernel_fwhm, size=kernel_size)
detconv = convolve(detimg, kernel, mask=coverage_mask)
WARNING: nan_treatment='interpolate', however, NaN values detected post convolution. A contiguous region of NaN values, larger than the kernel size, are present in the input array. Increase the kernel size to avoid this. [astropy.convolution.convolve]
coverage_mask = np.isnan(detimg)
bkg = Background2D(detimg, box_size=10, filter_size=5,
                   exclude_percentile=90, coverage_mask=coverage_mask, mask=ceers_bkgmask,
                   bkg_estimator=bkg_estimator, bkgrms_estimator=bkgrms_estimator)
bkgrms = bkg.background_rms
fig, ax = plt.subplots(ncols=2, figsize=(12, 6))
ax[0].imshow(detconv.value, norm=norm)
norm2 = simple_norm(bkgrms.value, 'sqrt', percent=99)
_ = ax[1].imshow(bkgrms.value, norm=norm2)
../../../_images/52fddaf595dd6d20a876a5048bc51190d8c1cc9d6b428379a47b34c015da2bca.png

Next, we need to define the detection threshold. In this example, we’ll define a 2D detection threshold image using the background RMS image calculated above. We set the threshold at the 4.0\(\sigma\) (per pixel) noise level (above the background).

threshold = bkgrms * 4.0
np.median(threshold)
\[0.0001850487 \; \mathrm{\mu Jy}\]

Now we are ready to detect the sources in the background-subtracted convolved image. Let’s find sources that have at least 25 connected pixels that are each greater than the corresponding threshold level defined above. For this, we use the detect_sources function.

mask = np.isnan(detimg)
npixels = 25
segm = detect_sources(detconv, threshold, npixels=npixels, mask=mask)
segm
<photutils.segmentation.core.SegmentationImage>
shape: (4800, 4700)
nlabels: 1785
labels: [   1    2    3    4    5 ... 1781 1782 1783 1784 1785]

The result is a segmentation image (SegmentationImage object). The segmentation image is an array with the same size as the science image, in which each detected source is labeled with a unique integer value (>= 1). Background pixels have a value of 0. As a simple example, a segmentation map containing two distinct sources (labeled 1 and 2) might look like this:

0 0 0 0 0 0 0 0 0 0
0 1 1 0 0 0 0 0 0 0
1 1 1 1 1 0 0 0 2 0
1 1 1 1 0 0 0 2 2 2
1 1 1 0 0 0 2 2 2 2
1 1 1 1 0 0 0 2 2 0
1 1 0 0 0 0 2 2 0 0
0 1 0 0 0 0 2 0 0 0
0 0 0 0 0 0 0 0 0 0

where all the pixels labeled 1 belong to the first source, all those labeled 2 belong to the second, and all those labeled 0 are background pixels.

Let’s plot the segmentation image.

fig, ax = plt.subplots(figsize=(8, 8))
_ = segm.imshow(ax=ax)
../../../_images/b41efb112c32184346b671aeb9ba36ddf1cfa8d265943e14c3104331f74f8fdd.png

Let’s plot the segmentation image next to the science image.

# switch to interactive plots
%matplotlib widget
fig, ax = plt.subplots(ncols=2, figsize=(12, 6), sharex=True, sharey=True)
ax[0].imshow(detconv.value, norm=norm)
ax[0].set_title('Detection Image')
segm.imshow(ax=ax[1])
_ = ax[1].set_title('Segmentation Image')

We can use the plot_patches SegmentationImage method to plot the segmentation outlines on the science image.

This is a good way to inspect the source extraction. The source detection parameters (e.g., threshold, npixels) can be adjusted as needed.

fig, ax = plt.subplots(figsize=(8, 8))
ax.imshow(detimg.value, norm=norm)
_ = segm.plot_patches(ax=ax)

The SegmentationImage object has attributes that can give various properties of the source segments.

star_label = segm.data[3852, 612]
idx = segm.get_index(label=star_label)
sobj = segm.segments[idx]
sobj.area, sobj.bbox, sobj.slices
(np.int64(149684),
 BoundingBox(ixmin=106, ixmax=1251, iymin=3297, iymax=4438),
 (slice(3297, 4438, None), slice(106, 1251, None)))
sobj.polygon
../../../_images/d807dc95fbce905ac1d4a2e843424beeef7d973cdd3e8c91443010400d9907ce.svg

Source deblending#

Comparing the data array with the segmentation image, we see that several detected sources were blended together. We can deblend them using the deblend_sources function, which uses a combination of multi-thresholding and watershed segmentation.

The deblending can be controlled with the deblend_sources keywords:

  • nlevels is the number of multi-thresholding levels to use; larger values give more deblending, but longer runtime

  • contrast is the fraction (0 - 1) of the total source flux that a local peak must have to be considered as a separate object

    • contrast = 1: no deblending

    • contrast = 0: max deblending

  • mode is the mode used in defining the spacing of the nlevels (‘exponential’, ‘linear’, ‘sinh’)

The npixels and connectivity values should match those used in detect_sources.

npixels2 = 9
segm_deblended = deblend_sources(detconv, segm, npixels=npixels2, progress_bar=True)
segm.nlabels, segm_deblended.nlabels, segm_deblended
(1785,
 2172,
 <photutils.segmentation.core.SegmentationImage>
 shape: (4800, 4700)
 nlabels: 2172
 labels: [   1    2    3    4    5 ... 2168 2169 2170 2171 2172])
fig, ax = plt.subplots(ncols=2, figsize=(12, 6), sharex=True, sharey=True)
segm.imshow(ax=ax[0])
ax[0].set_title('Original Segmentation Image')
segm_deblended.imshow(ax=ax[1])
_ = ax[1].set_title('Deblended Segmentation Image')
label = segm.data[1597, 3877]
idx = segm.get_index(label=label)
slc2 = segm.slices[idx]

fig, ax = plt.subplots(ncols=2, nrows=2, figsize=(10, 6))
ax = ax.ravel()
norm3 = simple_norm(detconv[slc2].value, 'sqrt', percent=98)
ax[0].imshow(detconv[slc2].value, norm=norm3)
ax[0].set_title('Detection Image')
segm[slc2].imshow(ax=ax[1])
ax[1].set_title('Segmentation Image')
segm_deblended[slc2].imshow(ax=ax[2])
ax[2].set_title('Deblended Segmentation Image')
ax[3].imshow(detimg[slc2].value, norm=norm3)
_ = segm_deblended[slc2].imshow(ax=ax[3], alpha=0.2)
child_labels = segm_deblended.deblended_labels_inverse_map[label]
child_labels
array([1749, 1750, 1751], dtype=int32)
segm_deblended.deblended_labels_map[child_labels[0]]
np.int32(577)

Modifying a Segmentation Image#

The SegmentationImage object provides several methods that can be used to modify itself (e.g., combining labels, removing labels, removing border segments) prior to measuring source photometry and other source properties, including:

  • reassign_labels: Reassign one or more label numbers.

  • relabel_consecutive: Reassign the label numbers consecutively, such that there are no missing label numbers.

  • keep_labels: Keep only the specified labels.

  • remove_labels: Remove one or more labels.

  • remove_border_labels: Remove labeled segments near the image border.

  • remove_masked_labels: Remove labeled segments located within a masked region.

Let’s combine the deblended segments of the bright star.

star_label = segm.data[3852, 612]
child_labels = segm_deblended.deblended_labels_inverse_map[star_label]
child_labels
array([2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013,
       2014, 2015, 2016, 2017, 2018, 2019], dtype=int32)

But, we’ll exclude the labels which are galaxies.

xy_exclude = ((462, 4059), (495, 3770), (807, 4100), (1110, 4270))
x_exclude, y_exclude = np.transpose(xy_exclude)
exclude_labels = segm_deblended.data[y_exclude, x_exclude]
exclude_labels
array([2014, 2004, 2017, 2019], dtype=int32)
keep_labels = np.setdiff1d(child_labels, exclude_labels)
keep_labels
array([2003, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013, 2015,
       2016, 2018], dtype=int32)
segm_deblended.reassign_labels(keep_labels, keep_labels[0])
fig, ax = plt.subplots(ncols=2, figsize=(12, 6), sharex=True, sharey=True)
segm.imshow(ax=ax[0])
ax[0].set_title('Original Segmentation Image')
segm_deblended.imshow(ax=ax[1])
_ = ax[1].set_title('Deblended Segmentation Image')

Let’s use the remove_label method to remove the bright star from the segmentation image.

segm_deblended.remove_label(keep_labels[0])
_ = segm_deblended.imshow(figsize=(8, 8))

When we combined and then removed segments, our segment labels are no longer consecutive.

segm_deblended
<photutils.segmentation.core.SegmentationImage>
shape: (4800, 4700)
nlabels: 2159
labels: [   1    2    3    4    5 ... 2168 2169 2170 2171 2172]
segm_deblended.is_consecutive
np.False_
segm_deblended.missing_labels
array([2003, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013, 2015,
       2016, 2018])
segm_deblended.relabel_consecutive()
segm_deblended.is_consecutive
np.True_
segm_deblended
<photutils.segmentation.core.SegmentationImage>
shape: (4800, 4700)
nlabels: 2159
labels: [   1    2    3    4    5 ... 2155 2156 2157 2158 2159]

SourceFinder class#

The SourceFinder class is a convenience class that combines the functionality of detect_sources and deblend_sources. After defining the SourceFinder object with the desired detection and deblending parameters, you call it with the background-subtracted (convolved) image and threshold.

# finder = SourceFinder(npixels=(npixels, npixels2), deblend=True, progress_bar=True)
# segment_map = finder(detconv, threshold, mask=mask)

Source Catalog: Photometry, Centroids, and Shape Properties#

The SourceCatalog class is used to measure the photometry, centroids, and shape/morphological properties of sources defined in a segmentation image. In its most basic form, it takes as input the (background-subtracted) image and the segmentation image. Usually the convolved image is also input, from which the source centroids and shape/morphological properties are measured (if not input, the unconvolved image is used instead).

First, we create a catalog from the detection image.

det_cat = SourceCatalog(detimg, segm_deblended, convolved_data=detconv, wcs=f356w['wcs'])
det_cat
<photutils.segmentation.catalog.SourceCatalog>
Length: 2159
labels: [   1    2    3    4    5 ... 2155 2156 2157 2158 2159]

Our det_cat variable is a SourceCatalog object for the detection image.

The properties of each source can be accessed using SourceCatalog attributes, or they can be output to an Astropy QTable using the to_table method.

Please see the SourceCatalog documentation for a complete list of the available source properties.

The to_table method is used to create a table of source measurements. The output is an Astropy QTable.

Each row in the table represents a source. The columns represent the calculated source properties. The label column corresponds to the label value in the input segmentation image.

det_tbl = det_cat.to_table()
det_tbl
QTable length=2159
labelxcentroidycentroidsky_centroidbbox_xminbbox_xmaxbbox_yminbbox_ymaxareasemimajor_sigmasemiminor_sigmaorientationeccentricitymin_valuemax_valuelocal_backgroundsegment_fluxsegment_fluxerrkron_fluxkron_fluxerr
deg,degpix2pixpixdeguJyuJyuJyuJyuJyuJyuJy
int32float64float64SkyCoordint64int64int64int64float64float64float64float64float64float64float64float64float64float64float64float64
1519.594955147681.07259985155056214.9609852444238,52.93199215314523514526798570.02.87246926378202261.6659138838634198-11.320560368908790.81464579262755320.000191590710587533580.0014283665865647930.00.04308027658363421nan0.06576522534470572nan
2766.48910486434398.5324969870587214.958588880522,52.93051982292866675178384114616.03.2239520142899023.06673131167734257.5926541360254410.30847161622966529.359696794881929e-050.240566290937736630.02.95442525875817nan2.8678235574541024nan
31371.252717570984396.86878719989326214.95318743218436,52.926673183264361348139389113758.07.857702791224084.288892164291958-18.45469530201170.83790238855002120.00016725764130857660.0105489322952791470.01.1038019934971528nan1.1939369727423381nan
43991.236136233889391.39582056328483214.9297798430202,52.9100151238763239883995909431.01.73666709486890871.2098242412307103-0.81451255492586650.7174250807358560.000167570730255129640.00149333416233354120.00.02059214366044146nan0.029888712810003045nan
51845.837111378475494.04547140465078214.94896532719426,52.9236461876953718421850939626.01.88165454237998091.0085658346537048-3.94688218297308070.84421831876973420.00024173994187358640.0011502134389694660.00.01713613686714633nan0.026092613208891995nan
63635.188122487194599.55083655878093214.93288180058389,52.91231919315401363236389310561.03.08851603800157951.4286982806507114-84.773795563089780.88657542809704790.000174248334519306360.0011094780990358390.00.03192307777654567nan0.05332255562080552nan
7295.4795965787447102.4665003529108214.962768593605,52.933529791885829329810010528.01.62967150695762061.2709620049314665-44.634400992845090.62591851468529720.00026872654612406060.00048936519184583130.00.009962071053491717nan0.04345487771200407nan
84205.705778270739102.9133784138537214.92773856854785,52.908715906394134203420910010636.01.8981929219786751.3299138437998494-45.341601350982740.71353360057153090.000182949450827979530.00096967280637508560.00.0195639354282886nan0.039308957758529224nan
93307.19135160857103.9118439804293214.93577283474275,52.914424647219723304331010110737.01.69500168293977121.4208569341215362-45.79765234198780.54526635840096330.000190694963425423820.00130277489992480120.00.022784642264923285nan0.04245106751305349nan
............................................................
21511318.20104556396384316.597489502455214.9092191179457,52.94979345128564131513214314431928.01.96084342762172351.085512184100727642.739482604577020.8327863488996880.00017341288537116930.000271739653278333940.00.006278585485429527nan0.06936167100627669nan
21521326.96103371882834321.127609327313214.90909290798857,52.94976228281462132213324317432547.02.99192266159296861.276695637320907737.345595085403990.90438656677792640.000121421705328417040.000381474039380452430.00.011397181107809979nan0.11311127078505667nan
2153494.186294544549644356.882692399877214.91617701196415,52.9552425216697148350943484367322.04.8024337098346512.967285601212213-29.238887798973830.78627930712883118.164345073692908e-050.0048153542539901440.00.2230271656583786nan0.2522678067290374nan
2154511.984351638372234351.301582630822214.91607635443594,52.9550994008455555095154349435430.01.5317627420665211.4459883452397193-12.8160757441401540.32993738699187210.00016017383493287190.000447937478628947760.00.007743445295495495nan0.02312506126273308nan
2155157.785609450891974423.995466755786214.9184840929336,52.95774047176380613417344064431678.06.0056260364786824.84632869559215-4.6716518055504060.59059956249618920.00012920178320833560.016609794791133120.01.4775174649673248nan1.6150349646531366nan
2156171.630131814690574424.209996831193214.91835777129148,52.9576537377104916618944074431447.05.01593989590726254.187121111769937-39.168851515448410.55060939729312390.000120931547201162090.0099002297590293030.00.7222632382815204nan0.8029248567175558nan
21574034.33319259576364467.067753876071214.88331213554403,52.93335746813930640084061444244711130.06.7765493839269355.161440263081692-2.68800486922331670.64797451638022235.105626814252766e-050.15956139803793080.01.934904479369893nan1.8738664801357052nan
21584002.5806512097524454.295761364674214.88373092850094,52.933490227483496400040054453445617.01.36479068349564050.93262065866243675.71137467010177250.73009714766456360.000109220462538143110.00036179636166117410.00.004565197407073732nan0.036476416958664nan
21594009.61547658057864453.691872786995214.88367432364515,52.93344229031217400640144451445633.02.06779504147073561.3282068385967711-2.9327341985891630.76642817033266150.000104255886024768120.00036877160226658880.00.008087272713176982nan0.066641109156109nan

The output columns can be customized using the optional columns keyword.

columns = ('label', 'xcentroid', 'ycentroid', 'area', 'segment_flux', 'kron_flux')
det_tbl2 = det_cat.to_table(columns=columns)
det_tbl2
QTable length=2159
labelxcentroidycentroidareasegment_fluxkron_flux
pix2uJyuJy
int32float64float64float64float64float64
1519.594955147681.0725998515505670.00.043080276583634210.06576522534470572
2766.48910486434398.5324969870587616.02.954425258758172.8678235574541024
31371.252717570984396.86878719989326758.01.10380199349715281.1939369727423381
43991.236136233889391.3958205632848331.00.020592143660441460.029888712810003045
51845.837111378475494.0454714046507826.00.017136136867146330.026092613208891995
63635.188122487194599.5508365587809361.00.031923077776545670.05332255562080552
7295.4795965787447102.466500352910828.00.0099620710534917170.04345487771200407
84205.705778270739102.913378413853736.00.01956393542828860.039308957758529224
93307.19135160857103.911843980429337.00.0227846422649232850.04245106751305349
..................
21511318.20104556396384316.59748950245528.00.0062785854854295270.06936167100627669
21521326.96103371882834321.12760932731347.00.0113971811078099790.11311127078505667
2153494.186294544549644356.882692399877322.00.22302716565837860.2522678067290374
2154511.984351638372234351.30158263082230.00.0077434452954954950.02312506126273308
2155157.785609450891974423.995466755786678.01.47751746496732481.6150349646531366
2156171.630131814690574424.209996831193447.00.72226323828152040.8029248567175558
21574034.33319259576364467.0677538760711130.01.9349044793698931.8738664801357052
21584002.5806512097524454.29576136467417.00.0045651974070737320.036476416958664
21594009.61547658057864453.69187278699533.00.0080872727131769820.066641109156109

Here, we access a few of the source properties as attributes of the SourceCatalog object. These get returned as NumPy arrays or Astropy objects (e.g., SkyCoord, Quantity).

det_cat.xcentroid, det_cat.sky_centroid, det_cat.segment_flux, det_cat.orientation
(array([ 519.59495515,  766.48910486, 1371.25271757, ..., 4034.3331926 ,
        4002.58065121, 4009.61547658], shape=(2159,)),
 <SkyCoord (ICRS): (ra, dec) in deg
     [(214.96098524, 52.93199215), (214.95858888, 52.93051982),
      (214.95318743, 52.92667318), ..., (214.88331214, 52.93335747),
      (214.88373093, 52.93349023), (214.88367432, 52.93344229)]>,
 <Quantity [0.04308028, 2.95442526, 1.10380199, ..., 1.93490448, 0.0045652 ,
            0.00808727] uJy>,
 <Quantity [-11.32056037,   7.59265414, -18.4546953 , ...,  -2.68800487,
              5.71137467,  -2.9327342 ] deg>)

Multiband Catalog#

Now let’s perform the photometry on each of the filter images. We pass in the science and (total) error array for each band, plus the segmentation image and the detection catalog SourceCatalog object created above. We need to input the detection catalog to ensure that the aperture centroids (and shapes for Kron apertures) are the same for all the bands.

f115w_cat = SourceCatalog(f115w['sci_bksub'], segm_deblended, error=f115w['err'], detection_cat=det_cat)
f200w_cat = SourceCatalog(f200w['sci_bksub'], segm_deblended, error=f200w['err'], detection_cat=det_cat)
f356w_cat = SourceCatalog(f356w['sci_bksub'], segm_deblended, error=f356w['err'], detection_cat=det_cat)

By default, the SourceCatalog class computes only “segment” (isophotal) photometry and Kron (elliptical aperture) photometry. We can also compute circular aperture photometry at various radii using the circular_photometry method.

f115w_cat.circular_photometry(radius=5, name='circ_rad5')
f115w_cat.circular_photometry(radius=10, name='circ_rad10')
f200w_cat.circular_photometry(radius=5, name='circ_rad5')
f200w_cat.circular_photometry(radius=10, name='circ_rad10')
f356w_cat.circular_photometry(radius=5, name='circ_rad5')
_ = f356w_cat.circular_photometry(radius=10, name='circ_rad10')
columns = ['segment_flux', 'segment_fluxerr', 'kron_flux', 'kron_fluxerr', 
           'circ_rad5_flux', 'circ_rad5_fluxerr', 'circ_rad10_flux', 'circ_rad10_fluxerr']
f115w_tbl = f115w_cat.to_table(columns=columns)
f200w_tbl = f200w_cat.to_table(columns=columns)
f356w_tbl = f356w_cat.to_table(columns=columns)
def add_magnitudes(tbl):
    for col in tbl.colnames:
        if col.endswith('flux'):
            fluxerr_col = col.replace('flux', 'fluxerr')
            flux = tbl[col]
            fluxerr = tbl[fluxerr_col]

            mag = -2.5 * np.log10(flux.value) + 23.9
            magerr = -2.5 * np.log10(1.0 + (flux / fluxerr))

            # set negative flux to the 2-sigma upper limit
            idx = flux.value <= 0
            mag[idx] = -2.5 * np.log10(2.0 * fluxerr[idx].value) + 23.9
            magerr[idx] = np.nan

            magcol = col.replace('flux', 'mag')
            magerrcol = f'{magcol}err'
            tbl[magcol] = mag * u.mag
            tbl[magerrcol] = magerr * u.mag
    return tbl


f115w_tbl = add_magnitudes(f115w_tbl)
f200w_tbl = add_magnitudes(f200w_tbl)
f356w_tbl = add_magnitudes(f356w_tbl)
/tmp/ipykernel_2196/270370035.py:8: RuntimeWarning: invalid value encountered in log10
  mag = -2.5 * np.log10(flux.value) + 23.9
/opt/hostedtoolcache/Python/3.11.13/x64/lib/python3.11/site-packages/astropy/units/quantity.py:659: RuntimeWarning: invalid value encountered in log10
  result = super().__array_ufunc__(function, method, *arrays, **kwargs)
def prefix_flux_columns(tbl, prefix):
    for col in tbl.colnames:
        if col != 'label':
            tbl.rename_column(col, f'{prefix}_{col}')
    return tbl


f115w_tbl = prefix_flux_columns(f115w_tbl, 'f115w')
f200w_tbl = prefix_flux_columns(f200w_tbl, 'f200w')
f356w_tbl = prefix_flux_columns(f356w_tbl, 'f356w')
f115w_tbl
QTable length=2159
f115w_segment_fluxf115w_segment_fluxerrf115w_kron_fluxf115w_kron_fluxerrf115w_circ_rad5_fluxf115w_circ_rad5_fluxerrf115w_circ_rad10_fluxf115w_circ_rad10_fluxerrf115w_segment_magf115w_segment_magerrf115w_kron_magf115w_kron_magerrf115w_circ_rad5_magf115w_circ_rad5_magerrf115w_circ_rad10_magf115w_circ_rad10_magerr
uJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmag
float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64
0.0264604267405598030.00327069847844541070.123249225518044930.0090125833776380870.0421915185848338660.00360350162898989950.120367231125483950.00732937938986090927.843507890018323-2.3964270928254426.173039502729186-2.91646290645328727.336937107311787-2.760233830823650726.1987293241901-3.1027804180285443
6.2446905405842770.0125325741246342666.0670222000805070.0113314180058166195.5027659513993170.0093933887987095735.9905894184437950.01094611422938101621.911222696205645-6.74585339750552921.94256103977239-6.82375425422438622.048547396500815-6.92124866199333821.956326112397267-6.847506014417105
1.20718976292163680.0103857340291142461.5651550568620510.0166268963502571450.42649264334062190.00367908650114616130.94226095280154820.00717234333896413723.695561140183333-5.17264675782045823.413606578195072-4.945838567098677524.825221139166867-5.16975466747241723.96457201407072-5.304508467844126
0.013826325137306290.00238982588052749630.024908944077938710.005646022637374090.0206686024316925580.00382190545559435740.021502339549122730.00756291416184657628.54823308636185-2.078953800523197327.90911170593235-1.83334799967924228.111722221194732-2.016794630456830328.068785710909054-1.4617123675908354
0.0193783080422890720.0022509479895234110.051164640113555570.0054819386416759880.0448197333033539540.003887223371272220.0564844955255583840.00755061312365926128.181710361811596-2.45669030310551227.12757519076435-2.535598787377465727.271326830079495-2.744878743731320527.020176864157296-2.321089830406173
0.031656386150297990.0029623522423207760.042665653911819010.0081385554643243280.0310521094639032440.00325733119225489820.047039799911604880.00704789531539265527.64884666276858-2.66918634818319627.32480398449952-1.98838092409462227.669772228862257-2.556379299332710427.218836335605236-2.212597577589724
0.0112258597011486750.00164542766287922860.062803362553510390.0062293927124158070.0237612387933160.00272625326179872860.0564092895626735340.00549839154783305628.774450974549215-2.23335798021755926.9050427577215-2.611523742562698527.96032730284277-2.468686603146154727.021623424722748-2.628772183884135
0.0231466469754837740.0020956660155206920.0305840784159067130.0056161019469258670.0269264706954277060.00304857686705601980.0324204737118583140.00607137114193793427.988779780028956-2.702017691899291327.68626150370067-2.02318937895755227.824551416689108-2.48165682762745527.62295160936473-2.005204861528172
0.0292566128044978830.0021194205619394780.044546971860383960.0059267973030945830.032467352434820990.00302718030944933770.044597727412829910.00610553990127971127.734439889773824-2.925952236907438727.277954530920887-2.325614156486473427.621382811683215-2.672807933999617427.276718178175585-2.298279680798929
................................................
0.0080030891174067180.00147111352998763320.106869739821205540.0066066855558520290.0180132807119688450.0024672249558948780.057973784282401660.00490526333670554429.14185586761688-2.022241200078215526.327863119790294-3.087305023743698828.26101795763054-2.29782981497208726.99192087451193-2.7696090698051017
0.0137126397964959290.00190441904123872520.160504660012044660.0085396182207760890.020118167173818190.00243980809983392580.055277785282099990.00493247217629042128.557197329712213-2.28459182855675325.886280884968066-3.241405050668882728.141028968403386-2.41486111876828827.043623413498125-2.7165146003532454
0.18558142236445020.0049776355735957620.184604369628622260.00735662106763076450.139191351323596640.0026280187555453150.172031339140558630.00492065285726937325.728663752518436-3.957516248979200425.7343950583411-3.54133649023383626.040969372106876-4.33026772919073825.81098107494007-3.8895818212812845
0.0091322956775240160.00147084088530391450.0098860215109448920.00524115075612090540.0104372572360780080.00237854054709000670.0109660570673154430.00477784111330374128.998550089363206-2.144671653448344528.91244612223754-1.150827751842412428.853534032006266-1.828587717319930328.79987374594994-1.2947014229665128
1.64063018387616770.0088083827868103983.21165570464769520.0133094707292206840.65472060959643750.00373860250556053851.67085244453961530.00681766021512644423.362473256694003-5.68109991332698822.633177545147518-5.96091060958422524.359859969853517-5.61454896324650323.34265475400941-5.977678032567718
0.71756291811158990.0070142010226845741.89999857007632020.0113274451774364480.335132811357943260.0033716726070168770.78498046155047070.00623207663394065924.26035003163938-5.035266018314783523.203116814733974-5.56800700527606525.08695762562139-5.00429757087766324.162852882170043-5.259150939340818
1.8439815764185870.00971963349729776414.02311727295670.0173446730470575332.29810769784195920.00582432744305408212.5896648101115540.01449421875580428223.235608555955043-5.7009746107670321.032888585046216-7.27053819779408122.996574056305484-6.49305969657354521.14996458121501-7.348297656075529
0.0057046917366317590.00104165694210678340.0491741527394816340.0045761492815864180.0161163300421351060.00225810390685679630.046092083943435240.00460296380239325629.509419546255028-2.028360177466824527.170657785371475-2.67470181364246228.38183461847187-2.27617511685980827.24093413961089-2.604819941987289
0.0082551306243779010.00143856811337172990.073083387062537970.0064590781198572920.0134998965072163150.0021583860513611790.033385724894412410.00421312217386239329.108190125809845-2.07139772205802226.740453333495687-2.72607127335031628.574173902191585-2.151537500208593527.591097974049852-2.3764261868930676

Now let’s combine the detection and filter tables into one main table.

det_columns = det_tbl.colnames[:-7]  # trim some columns
main_tbl = det_cat.to_table(det_columns)
main_tbl
QTable length=2159
labelxcentroidycentroidsky_centroidbbox_xminbbox_xmaxbbox_yminbbox_ymaxareasemimajor_sigmasemiminor_sigmaorientationeccentricity
deg,degpix2pixpixdeg
int32float64float64SkyCoordint64int64int64int64float64float64float64float64float64
1519.594955147681.07259985155056214.9609852444238,52.93199215314523514526798570.02.87246926378202261.6659138838634198-11.320560368908790.8146457926275532
2766.48910486434398.5324969870587214.958588880522,52.93051982292866675178384114616.03.2239520142899023.06673131167734257.5926541360254410.3084716162296652
31371.252717570984396.86878719989326214.95318743218436,52.926673183264361348139389113758.07.857702791224084.288892164291958-18.45469530201170.8379023885500212
43991.236136233889391.39582056328483214.9297798430202,52.9100151238763239883995909431.01.73666709486890871.2098242412307103-0.81451255492586650.717425080735856
51845.837111378475494.04547140465078214.94896532719426,52.9236461876953718421850939626.01.88165454237998091.0085658346537048-3.94688218297308070.8442183187697342
63635.188122487194599.55083655878093214.93288180058389,52.91231919315401363236389310561.03.08851603800157951.4286982806507114-84.773795563089780.8865754280970479
7295.4795965787447102.4665003529108214.962768593605,52.933529791885829329810010528.01.62967150695762061.2709620049314665-44.634400992845090.6259185146852972
84205.705778270739102.9133784138537214.92773856854785,52.908715906394134203420910010636.01.8981929219786751.3299138437998494-45.341601350982740.7135336005715309
93307.19135160857103.9118439804293214.93577283474275,52.914424647219723304331010110737.01.69500168293977121.4208569341215362-45.79765234198780.5452663584009633
.......................................
21511318.20104556396384316.597489502455214.9092191179457,52.94979345128564131513214314431928.01.96084342762172351.085512184100727642.739482604577020.832786348899688
21521326.96103371882834321.127609327313214.90909290798857,52.94976228281462132213324317432547.02.99192266159296861.276695637320907737.345595085403990.9043865667779264
2153494.186294544549644356.882692399877214.91617701196415,52.9552425216697148350943484367322.04.8024337098346512.967285601212213-29.238887798973830.7862793071288311
2154511.984351638372234351.301582630822214.91607635443594,52.9550994008455555095154349435430.01.5317627420665211.4459883452397193-12.8160757441401540.3299373869918721
2155157.785609450891974423.995466755786214.9184840929336,52.95774047176380613417344064431678.06.0056260364786824.84632869559215-4.6716518055504060.5905995624961892
2156171.630131814690574424.209996831193214.91835777129148,52.9576537377104916618944074431447.05.01593989590726254.187121111769937-39.168851515448410.5506093972931239
21574034.33319259576364467.067753876071214.88331213554403,52.93335746813930640084061444244711130.06.7765493839269355.161440263081692-2.68800486922331670.6479745163802223
21584002.5806512097524454.295761364674214.88373092850094,52.933490227483496400040054453445617.01.36479068349564050.93262065866243675.71137467010177250.7300971476645636
21594009.61547658057864453.691872786995214.88367432364515,52.93344229031217400640144451445633.02.06779504147073561.3282068385967711-2.9327341985891630.7664281703326615
mult_tbl = hstack((main_tbl, f115w_tbl, f200w_tbl, f356w_tbl), metadata_conflicts='silent')
mult_tbl
QTable length=2159
labelxcentroidycentroidsky_centroidbbox_xminbbox_xmaxbbox_yminbbox_ymaxareasemimajor_sigmasemiminor_sigmaorientationeccentricityf115w_segment_fluxf115w_segment_fluxerrf115w_kron_fluxf115w_kron_fluxerrf115w_circ_rad5_fluxf115w_circ_rad5_fluxerrf115w_circ_rad10_fluxf115w_circ_rad10_fluxerrf115w_segment_magf115w_segment_magerrf115w_kron_magf115w_kron_magerrf115w_circ_rad5_magf115w_circ_rad5_magerrf115w_circ_rad10_magf115w_circ_rad10_magerrf200w_segment_fluxf200w_segment_fluxerrf200w_kron_fluxf200w_kron_fluxerrf200w_circ_rad5_fluxf200w_circ_rad5_fluxerrf200w_circ_rad10_fluxf200w_circ_rad10_fluxerrf200w_segment_magf200w_segment_magerrf200w_kron_magf200w_kron_magerrf200w_circ_rad5_magf200w_circ_rad5_magerrf200w_circ_rad10_magf200w_circ_rad10_magerrf356w_segment_fluxf356w_segment_fluxerrf356w_kron_fluxf356w_kron_fluxerrf356w_circ_rad5_fluxf356w_circ_rad5_fluxerrf356w_circ_rad10_fluxf356w_circ_rad10_fluxerrf356w_segment_magf356w_segment_magerrf356w_kron_magf356w_kron_magerrf356w_circ_rad5_magf356w_circ_rad5_magerrf356w_circ_rad10_magf356w_circ_rad10_magerr
deg,degpix2pixpixdeguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmaguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmaguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmag
int32float64float64SkyCoordint64int64int64int64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64
1519.594955147681.07259985155056214.9609852444238,52.93199215314523514526798570.02.87246926378202261.6659138838634198-11.320560368908790.81464579262755320.0264604267405598030.00327069847844541070.123249225518044930.0090125833776380870.0421915185848338660.00360350162898989950.120367231125483950.00732937938986090927.843507890018323-2.3964270928254426.173039502729186-2.91646290645328727.336937107311787-2.760233830823650726.1987293241901-3.10278041802854430.043236918493417520.00329770613461732860.07997834108295451nan0.040387477885375854nan0.06604594846733264nan27.310363164745095-2.87391065716325726.642569020698495nan27.384383166942825nan26.85038454664706nan0.048500640877349510.00204403907991945740.07874897788427716nan0.055891662886592156nan0.07781177288312524nan27.18563130670787-3.482965640462023526.659387685976032nan27.031632422698127nan26.672386723741447nan
2766.48910486434398.5324969870587214.958588880522,52.93051982292866675178384114616.03.2239520142899023.06673131167734257.5926541360254410.30847161622966526.2446905405842770.0125325741246342666.0670222000805070.0113314180058166195.5027659513993170.0093933887987095735.9905894184437950.01094611422938101621.911222696205645-6.74585339750552921.94256103977239-6.82375425422438622.048547396500815-6.92124866199333821.956326112397267-6.8475060144171055.1855555644039210.0170160867273807535.07956046126657950.0154969723766715984.6194839838294880.013216535341674095.03219187248708360.01507941494730063322.113011769532775-6.21339594932780622.135434664912825-6.29225556252882722.23851633572635-6.3617915113200622.145607020072646-6.3116803908467341.8374013717025010.0057437852956354621.78431590454676630.0048716169078053921.41684252104437540.0031165856107116271.76761917470223030.00458866541354715223.23948990935086-6.26590335492031323.271320633509827-6.41245682865304923.52169604488616-6.64649194034178123.281528239908603-6.467070654037891
31371.252717570984396.86878719989326214.95318743218436,52.926673183264361348139389113758.07.857702791224084.288892164291958-18.45469530201170.83790238855002121.20718976292163680.0103857340291142461.5651550568620510.0166268963502571450.42649264334062190.00367908650114616130.94226095280154820.00717234333896413723.695561140183333-5.17264675782045823.413606578195072-4.945838567098677524.825221139166867-5.16975466747241723.96457201407072-5.3045084678441261.26659884420245520.0125124063342809681.3544375424123594nan0.446305004381272460.0048583991704462330.9676378324972978nan23.643402281527955-5.02391869845582723.57060254230436nan24.775920608559737-4.9196016758256223.935717900063292nan0.99147753510519170.0069801392965018751.145851220507556nan0.31326643390398150.00275381279331680160.7458282143811255nan23.909292803975443-5.38866392471486823.75217942081515nan25.160215341803365-5.149451236437454524.218402978696297nan
43991.236136233889391.39582056328483214.9297798430202,52.9100151238763239883995909431.01.73666709486890871.2098242412307103-0.81451255492586650.7174250807358560.013826325137306290.00238982588052749630.024908944077938710.005646022637374090.0206686024316925580.00382190545559435740.021502339549122730.00756291416184657628.54823308636185-2.078953800523197327.90911170593235-1.83334799967924228.111722221194732-2.016794630456830328.068785710909054-1.46171236759083540.026953876152693090.00232771900482475760.035462554177720124nan0.031804712236071064nan0.03729965754914847nan27.823446928634834-2.74916044856232927.525574969261505nan27.643771323741007nan27.47073788865508nan0.020907546315298650.0021078586578369140.03838842741081264nan0.03171930952004332nan0.043248392034779176nan28.09924233118558-2.595442837672558227.439499196184418nan27.646690687914685nan27.31007508716868nan
51845.837111378475494.04547140465078214.94896532719426,52.9236461876953718421850939626.01.88165454237998091.0085658346537048-3.94688218297308070.84421831876973420.0193783080422890720.0022509479895234110.051164640113555570.0054819386416759880.0448197333033539540.003887223371272220.0564844955255583840.00755061312365926128.181710361811596-2.45669030310551227.12757519076435-2.535598787377465727.271326830079495-2.744878743731320527.020176864157296-2.3210898304061730.015616735137302920.00258455099537968640.024647970745112286nan0.021073677088702166nan0.029767593468605243nan28.416024388384116-2.119292428606185427.92054707539663nan28.090649197287068nan27.715640685095558nan0.0177843312248679440.00135283626150339840.03336803632887296nan0.03060621172791021nan0.03726895417919314nan28.274906154254325-2.87658105545807927.591673375828385nan27.68547605427185nan27.471631985017606nan
63635.188122487194599.55083655878093214.93288180058389,52.91231919315401363236389310561.03.08851603800157951.4286982806507114-84.773795563089780.88657542809704790.031656386150297990.0029623522423207760.042665653911819010.0081385554643243280.0310521094639032440.00325733119225489820.047039799911604880.00704789531539265527.64884666276858-2.66918634818319627.32480398449952-1.98838092409462227.669772228862257-2.556379299332710427.218836335605236-2.2125975775897240.039675927763107270.0035077773500233890.05892264085199979nan0.04047530505520450.0039305214291506340.06257033955498775nan27.403682272794654-2.725719705199673526.974294491905955nan27.38202467340418-2.632474473014820426.909078720497007nan0.029097616833617040.00189000566024333240.060733574303710416nan0.031626316429826010.0021587550362111610.06465084500371643nan27.740356448359353-3.036812879323984626.94142789736241nan27.649878470365582-2.98630369544633426.87356448600564nan
7295.4795965787447102.4665003529108214.962768593605,52.933529791885829329810010528.01.62967150695762061.2709620049314665-44.634400992845090.62591851468529720.0112258597011486750.00164542766287922860.062803362553510390.0062293927124158070.0237612387933160.00272625326179872860.0564092895626735340.00549839154783305628.774450974549215-2.23335798021755926.9050427577215-2.611523742562698527.96032730284277-2.468686603146154727.021623424722748-2.6287721838841350.0107178403640096460.0021007398609071970.06472108106714805nan0.0273840457414507160.00353941022527190850.053967083758949166inf28.824731789519042-1.963669124082960226.87238559276728nan27.806255971094245-2.3533928073629927.069677623616613-0.00.0093866398230337460.00098316383082419630.033939036830723190.00324930225913691320.0168589744633006850.0015595894714486220.032275174021595310.002944985143896961428.96872461557558-2.557861600624275627.573251217143493-2.646541674434347528.332922117943028-2.680613678204571727.627828518555106-2.6942705149244883
84205.705778270739102.9133784138537214.92773856854785,52.908715906394134203420910010636.01.8981929219786751.3299138437998494-45.341601350982740.71353360057153090.0231466469754837740.0020956660155206920.0305840784159067130.0056161019469258670.0269264706954277060.00304857686705601980.0324204737118583140.00607137114193793427.988779780028956-2.702017691899291327.68626150370067-2.02318937895755227.824551416689108-2.48165682762745527.62295160936473-2.0052048615281720.025658217491788370.00252682273276150230.0479529765908209560.0070267612078081420.034624681444209440.0036663188058932630.0509930073691267150.00753350815069754927.876933794724398-2.618609708898631527.197961074108886-2.233618673222047827.551535533819347-2.547166211997705427.131223435697017-2.22588851477110960.019212683642199902inf0.04457946470439106inf0.030443840224226942inf0.043220212295755386inf28.19102992079778-0.027.277162876725864-0.027.691251414987057-0.027.310782760531026-0.0
93307.19135160857103.9118439804293214.93577283474275,52.914424647219723304331010110737.01.69500168293977121.4208569341215362-45.79765234198780.54526635840096330.0292566128044978830.0021194205619394780.044546971860383960.0059267973030945830.032467352434820990.00302718030944933770.044597727412829910.00610553990127971127.734439889773824-2.925952236907438727.277954530920887-2.325614156486473427.621382811683215-2.672807933999617427.276718178175585-2.2982796807989290.0287202170965471030.00257594091817736630.048627387199156840.0074881743599930210.03214413700719690.00385505922276621660.047501388841803830.0076888136528994727.754530703947243-2.71138782013753827.182797662264676-2.186768412683106627.63224557377189-2.425654374920215527.208234230805523-2.14000663024004820.017354869846940730.00146078120451420550.0373965256373695250.0040249535670844380.0234298870241751380.00213008296821642860.039520387983429110.00412594750821851728.301446597442872-2.774835194804805627.467921861250748-2.53116182336866427.975574513771743-2.697909549451207827.407947001162743-2.561060220407469
.......................................................................................................................................................................................
21511318.20104556396384316.597489502455214.9092191179457,52.94979345128564131513214314431928.01.96084342762172351.085512184100727642.739482604577020.8327863488996880.0080030891174067180.00147111352998763320.106869739821205540.0066066855558520290.0180132807119688450.0024672249558948780.057973784282401660.00490526333670554429.14185586761688-2.022241200078215526.327863119790294-3.087305023743698828.26101795763054-2.29782981497208726.99192087451193-2.76960906980510170.0058589186024572490.00121425068937242030.100030104144745590.0055811215620171850.0147322754517803090.0020497660205130460.0555741274139381350.00412837147164599729.48045633893291-1.913264235152765126.39967319758943-3.192471502613277528.479325424164404-2.282851251380153527.037818368908763-2.9005343522526470.0061597606333856890.00055167858954519030.057271582371854080.0024736199419084060.0124674302324843420.00091887332961776930.027671735111603340.00182923524021817629.426090410181292-2.712823830539290427.005152043706513-3.457425156617783428.6605576335173-2.908512557579969327.79490902067226-3.018916854103474
21521326.96103371882834321.127609327313214.90909290798857,52.94976228281462132213324317432547.02.99192266159296861.276695637320907737.345595085403990.90438656677792640.0137126397964959290.00190441904123872520.160504660012044660.0085396182207760890.020118167173818190.00243980809983392580.055277785282099990.00493247217629042128.557197329712213-2.28459182855675325.886280884968066-3.241405050668882728.141028968403386-2.41486111876828827.043623413498125-2.71651460035324540.0110581545202680510.00167247024364769460.160532783789712780.007351289135269780.0173712554681190940.00214212644019360620.055797008261188660.0041643938429384428.79079336443621-2.20372829466365625.886090657913268-3.39661498094147728.300421981923016-2.398718606880620627.033472716387656-2.8957999597379580.0111921107715732410.00075574935181066390.095248418942681770.003288028800125820.0135117080413112990.00097626842831919710.034566387836077930.001919578828551462428.77772000039061-2.997280865503313426.45285556144497-3.69165322229337128.57322436874236-2.928596217818017627.55336500390812-3.197299764110899
2153494.186294544549644356.882692399877214.91617701196415,52.9552425216697148350943484367322.04.8024337098346512.967285601212213-29.238887798973830.78627930712883110.18558142236445020.0049776355735957620.184604369628622260.00735662106763076450.139191351323596640.0026280187555453150.172031339140558630.00492065285726937325.728663752518436-3.957516248979200425.7343950583411-3.54133649023383626.040969372106876-4.33026772919073825.81098107494007-3.88958182128128450.212508667532091130.0043986290693283080.23775248847812820.00650107459502042950.145721328885928180.00229878871421485820.206568036250614260.00432468573887623825.581558379509563-4.232392053935015525.459687321868447-3.937139396187719225.991192192305846-4.52205422148661325.612342197771024-4.2202675953299080.229299970055986070.00185626780148595570.263026111709060.0026950845560724250.14114841569570250.00099873054568954930.220625737340483960.00183711830587583925.49899000496016-5.238162465079958525.350002837822696-4.98463450894521526.02580998044649-5.38322453103533825.540859564588175-5.207800957623072
2154511.984351638372234351.301582630822214.91607635443594,52.9550994008455555095154349435430.01.5317627420665211.4459883452397193-12.8160757441401540.32993738699187210.0091322956775240160.00147084088530391450.0098860215109448920.00524115075612090540.0104372572360780080.00237854054709000670.0109660570673154430.00477784111330374128.998550089363206-2.144671653448344528.91244612223754-1.150827751842412428.853534032006266-1.828587717319930328.79987374594994-1.29470142296651280.0086866194868396240.0013286682078614830.0127141767690226920.004606158730095950.011431548914768060.0021629829792194520.013324446340112850.00422393062464460429.0528730053246-2.193117214379937328.639279386842034-1.43804350075713428.754737302497794-1.995777890019894128.5883770695301-1.5463004463426080.00730960671589063640.00054405181435868140.026473030794909870.00190039598624970220.0135077023830775920.00087142686087456460.025347533476558440.001729806585749311229.240264972691975-2.8985793791457427.842990837891733-2.93516921522389528.57354629199174-3.043754111982966727.890160736747383-2.9865211325968866
2155157.785609450891974423.995466755786214.9184840929336,52.95774047176380613417344064431678.06.0056260364786824.84632869559215-4.6716518055504060.59059956249618921.64063018387616770.0088083827868103983.21165570464769520.0133094707292206840.65472060959643750.00373860250556053851.67085244453961530.00681766021512644423.362473256694003-5.68109991332698822.633177545147518-5.96091060958422524.359859969853517-5.61454896324650323.34265475400941-5.9776780325677181.92688206044578750.0098355859518051153.701369648753726nan0.77504853620507370.0038746787324177911.98333876630720820.00704627348222567523.18786216648163-5.735665179170295522.479093851435763nan24.1766777490697-5.75814737182341323.15650774830886-6.1274440090316221.08869894588080870.0061898655258119111.1906367377131581nan0.402457566342197370.00284697194290164980.9282962262407563nan23.807730493949496-5.619222009114551523.71055180299092nan24.888199758270865-5.38349569873475623.98078353765037nan
2156171.630131814690574424.209996831193214.91835777129148,52.9576537377104916618944074431447.05.01593989590726254.187121111769937-39.168851515448410.55060939729312390.71756291811158990.0070142010226845741.89999857007632020.0113274451774364480.335132811357943260.0033716726070168770.78498046155047070.00623207663394065924.26035003163938-5.035266018314783523.203116814733974-5.56800700527606525.08695762562139-5.00429757087766324.162852882170043-5.2591509393408180.88696802654542310.00687686493620276452.28315180110905920.0110803930687160220.417658672662045370.0033434338264903670.96998829714322960.00618880734357675724.03023008837211-5.284679169137036523.003663030963065-5.79020550708971624.84794624015169-5.25022884202078523.933083763604433-5.4948041349672070.4930229402197040.0060713244602084160.5547875147067951nan0.21751037725613550.0041359385670294010.45834210128135155nan24.667832181629215-4.78724784102769924.53968330282419nan25.556300045881777-4.32271609386703324.747025621374252nan
21574034.33319259576364467.067753876071214.88331213554403,52.93335746813930640084061444244711130.06.7765493839269355.161440263081692-2.68800486922331670.64797451638022231.8439815764185870.00971963349729776414.02311727295670.0173446730470575332.29810769784195920.00582432744305408212.5896648101115540.01449421875580428223.235608555955043-5.7009746107670321.032888585046216-7.27053819779408122.996574056305484-6.49305969657354521.14996458121501-7.3482976560755292.23339992095671040.00909776054322719612.841343776271290.0156673083025696822.59450112335642120.006054852877867120511.6428825877983950.01333335985598269223.02758375892064-5.97949375216648321.128473818377515-7.28536407183917222.864865342324936-6.58240652813078321.23484870536059-7.3540449606047051.77450963597630060.0092895189300179481.717806938296575nan1.0960666189483448nan1.4893995368687927nan23.277304094989375-5.70838182443309923.312564118642513nan23.800407621557984nan23.467471963843007nan
21584002.5806512097524454.295761364674214.88373092850094,52.933490227483496400040054453445617.01.36479068349564050.93262065866243675.71137467010177250.73009714766456360.0057046917366317590.00104165694210678340.0491741527394816340.0045761492815864180.0161163300421351060.00225810390685679630.046092083943435240.00460296380239325629.509419546255028-2.028360177466824527.170657785371475-2.67470181364246228.38183461847187-2.27617511685980827.24093413961089-2.6048199419872890.0053093787160579940.00102793902624398470.040736364710168990.0042273752389510650.0156235708568706370.0021078331769543690.0416366301112147260.00422633689112783929.587390738715243-1.974845319664029227.375044324493086-2.56697907024788428.41554924653755-2.312267224863184427.351311069141033-2.58874508574420140.00392266772603787950.00065710110357031230.027160470351268960.00314571780648998560.0091571590067070440.0015043034772759040.0257623526229424570.003236360113924460829.916046194282202-2.108028400882054227.815156783586357-2.459528886377437428.99559811143271-2.126203309086079727.87253619911634-2.3808046992751866
21594009.61547658057864453.691872786995214.88367432364515,52.93344229031217400640144451445633.02.06779504147073561.3282068385967711-2.9327341985891630.76642817033266150.0082551306243779010.00143856811337172990.073083387062537970.0064590781198572920.0134998965072163150.0021583860513611790.033385724894412410.00421312217386239329.108190125809845-2.07139772205802226.740453333495687-2.72607127335031628.574173902191585-2.151537500208593527.591097974049852-2.37642618689306760.0057191084497626480.0013019035104662180.073005178668154790.0059233282749880630.0110009062151039550.00200141960634211270.0367529575720643850.00390776397039013129.50667916994346-1.829552287035034526.741615829588177-2.811670213271346328.79642884437525-2.03170723696045127.486769266922607-2.54311689184886670.0093224762881193390.00093827769160270690.057197350129722990.0045883786669067070.0164105695738174140.001448517196227080.0283354760239643220.00297030904549787628.976171781502053-2.597119710067710627.006560227435997-2.8230723047173528.362190863261052-2.727338979982683527.769173717495825-2.557057403665186
len(mult_tbl.colnames), mult_tbl.colnames
(61,
 ['label',
  'xcentroid',
  'ycentroid',
  'sky_centroid',
  'bbox_xmin',
  'bbox_xmax',
  'bbox_ymin',
  'bbox_ymax',
  'area',
  'semimajor_sigma',
  'semiminor_sigma',
  'orientation',
  'eccentricity',
  'f115w_segment_flux',
  'f115w_segment_fluxerr',
  'f115w_kron_flux',
  'f115w_kron_fluxerr',
  'f115w_circ_rad5_flux',
  'f115w_circ_rad5_fluxerr',
  'f115w_circ_rad10_flux',
  'f115w_circ_rad10_fluxerr',
  'f115w_segment_mag',
  'f115w_segment_magerr',
  'f115w_kron_mag',
  'f115w_kron_magerr',
  'f115w_circ_rad5_mag',
  'f115w_circ_rad5_magerr',
  'f115w_circ_rad10_mag',
  'f115w_circ_rad10_magerr',
  'f200w_segment_flux',
  'f200w_segment_fluxerr',
  'f200w_kron_flux',
  'f200w_kron_fluxerr',
  'f200w_circ_rad5_flux',
  'f200w_circ_rad5_fluxerr',
  'f200w_circ_rad10_flux',
  'f200w_circ_rad10_fluxerr',
  'f200w_segment_mag',
  'f200w_segment_magerr',
  'f200w_kron_mag',
  'f200w_kron_magerr',
  'f200w_circ_rad5_mag',
  'f200w_circ_rad5_magerr',
  'f200w_circ_rad10_mag',
  'f200w_circ_rad10_magerr',
  'f356w_segment_flux',
  'f356w_segment_fluxerr',
  'f356w_kron_flux',
  'f356w_kron_fluxerr',
  'f356w_circ_rad5_flux',
  'f356w_circ_rad5_fluxerr',
  'f356w_circ_rad10_flux',
  'f356w_circ_rad10_fluxerr',
  'f356w_segment_mag',
  'f356w_segment_magerr',
  'f356w_kron_mag',
  'f356w_kron_magerr',
  'f356w_circ_rad5_mag',
  'f356w_circ_rad5_magerr',
  'f356w_circ_rad10_mag',
  'f356w_circ_rad10_magerr'])

Let’s create a copy of our source table where we included the largest 500 sources.

tmp_tbl = mult_tbl.copy()
tmp_tbl.sort(keys='area', reverse=True)  # sort by decreasing area
tmp_tbl = tmp_tbl[0:500]  # select first 500 sources
tmp_tbl
QTable length=500
labelxcentroidycentroidsky_centroidbbox_xminbbox_xmaxbbox_yminbbox_ymaxareasemimajor_sigmasemiminor_sigmaorientationeccentricityf115w_segment_fluxf115w_segment_fluxerrf115w_kron_fluxf115w_kron_fluxerrf115w_circ_rad5_fluxf115w_circ_rad5_fluxerrf115w_circ_rad10_fluxf115w_circ_rad10_fluxerrf115w_segment_magf115w_segment_magerrf115w_kron_magf115w_kron_magerrf115w_circ_rad5_magf115w_circ_rad5_magerrf115w_circ_rad10_magf115w_circ_rad10_magerrf200w_segment_fluxf200w_segment_fluxerrf200w_kron_fluxf200w_kron_fluxerrf200w_circ_rad5_fluxf200w_circ_rad5_fluxerrf200w_circ_rad10_fluxf200w_circ_rad10_fluxerrf200w_segment_magf200w_segment_magerrf200w_kron_magf200w_kron_magerrf200w_circ_rad5_magf200w_circ_rad5_magerrf200w_circ_rad10_magf200w_circ_rad10_magerrf356w_segment_fluxf356w_segment_fluxerrf356w_kron_fluxf356w_kron_fluxerrf356w_circ_rad5_fluxf356w_circ_rad5_fluxerrf356w_circ_rad10_fluxf356w_circ_rad10_fluxerrf356w_segment_magf356w_segment_magerrf356w_kron_magf356w_kron_magerrf356w_circ_rad5_magf356w_circ_rad5_magerrf356w_circ_rad10_magf356w_circ_rad10_magerr
deg,degpix2pixpixdeguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmaguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmaguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmag
int32float64float64SkyCoordint64int64int64int64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64
19221845.44577999075272740.020735397627214.92110658578846,52.93793594738365176819272673281214527.022.30832783754874516.549664468454985-28.958220458709190.670554384910276514.897869719260020.0321815609931945814.9656565373052160.03183507352244750.90624252208214520.00327745574466014762.5312756193693360.00600316742720008420.967189569793668-6.666135483281420.962260565206073-6.68278192620581824.006888910036814-6.10818851686325922.891651409670374-6.56496933149338432.694880321458920.03207345679402351432.781326548626480.031825236094924562.62508914143890370.005137651773497716.9518278994168480.00863921451705019420.113800620047666-7.5218995337894220.110933689978804-7.533190783274875522.85213986114248-6.77307133178556421.794752469728003-7.26541031990858854.0819135008876940.01548276096582412754.2327979902222950.0153977180289635014.5964402079981620.003192338509888063512.5230160587495760.00528111184412906519.567369876977484-8.85831988471437419.56434497192965-8.8673223395334522.243945962388988-7.89653550350139521.155727656263352-8.437916702860873
17573683.74749806715451728.8248032182457214.91529410352754,52.920806831597034362437501639182912201.028.322546047159157.83253617627225561.405040319429430.961000180030673212.268927964145260.0338472202420234711.9184057136344810.0286500693942415421.05179478189250020.00371282825012529162.62536916500107330.00664142214842783121.177983458625896-6.40120020155985821.209454587177042-6.5503380331909823.845172470288894-6.13439127253014722.852024049644825-6.49506636429701736.923399678189470.0330621413886547136.2052425631218140.0291391069431195354.6953286702578610.00680971848283072210.3826190679645230.01045864918959074619.981745796640766-7.62089851689773220.00307134636156-7.73661154310419422.220835004935893-7.09791562093742721.359232699241705-7.49317144746649361.336479014193220.0162871628999710160.042591298459970.0149621135833489367.9911144274713070.00462701231840201218.402181806956360.00687943078760888319.430702895166345-8.93997177048462519.45385143221377-9.0089367221948921.64348162619828-8.09389521518328820.73782670734738-8.568697841831721
3203774.88416260942131115.8669131621625214.92093215787574,52.916919417029085371438371051118211366.014.86839995275105211.26726599653185859.273286532679620.652487488960095827.8612165491304540.03481952473521232626.830526130162810.0275702547562899088.7934607406768830.00847298339891911314.2618665542125150.01137254761150291620.287499810498076-7.25929914251198220.328427027490285-7.47158612229070221.539600427564398-7.54135434658882821.01455907871642-7.7466619433908753.336509648448790.0355520285665988951.5052644567310550.03039808188391728717.89500564622430.0150685036246020428.811232962176020.0183706230332541619.582438520078796-7.94112398250675419.62037094644267-8.0731542078562820.76817040061388-7.6875681424340920.25109538977174-7.98928696315625854.574279096950790.01653228700160980252.732484622603550.01471211977034814216.4419729110574480.00815877694899763929.1228965655707630.01023274891822677219.557529981066484-8.79696653305538819.594804413295428-8.8863103303309820.86011517912841-8.26136079685158720.23941358061906-8.63598704657102
17912619.7052472773961830.0070059551974214.92375635469517,52.928108055052576254526881765189111215.018.16438124932532716.1354427091301332.677177095203650.459261228025413523.9213752150686860.0314903222024440824.241017251953680.0315787210077567652.85805613204524130.0050060060567056675.8644462979837310.00791652156688162820.45303464206569-7.20295093062978920.438622898388918-7.21430422647208422.75982311486147-6.893348507645896521.979432465726248-7.17569620262300246.660523402024660.0321990065276622847.145694057670740.0323252452064701846.3810130568279780.00801570606969932713.1026747099118770.01159749781476922419.727625984829935-7.90351680996462519.716394915979862-7.910494689477709521.887775916808486-7.25373264938380321.106600101855864-7.63344973635447153.4499310370396050.01454388070851564453.9570674873437640.014592976219574266.78145001160009550.00400320441657469715.1926461537011440.00588106620218427319.580132126991725-8.91346250545008819.56987915457376-8.92005476198015221.821693588080468-8.07292773334489220.945916442587272-8.530863593537035
21404200.0408360278284342.8708913563005214.88313759904435,52.9316349823253412742674277439311033.013.16797797961782711.280399266028091-56.5838223469400250.515891607757506447.581280204304630.03098545782268047345.973572575422170.02511019420777134816.9096780214855720.01112873070343783927.088604381886670.01404101696485343719.706409692937633-7.96640232688859819.743729365515158-8.15723830908409320.82966165439763-7.954938580488799520.31803342363512-8.21403279725434277.860457271914830.03710841760039329575.297948899420490.0327043763545711128.9182014505505940.02003586011857025746.1019796786547150.02360583625195418819.171707625032756-8.30512862705548119.208042134325463-8.40591465447176820.247071803028035-7.8991602070349219.740701062607453-8.2273062547596642.656546544225270.01436535920947790141.1499927559181060.01228136230539663213.8515797550266250.006565498617242853524.3737899762360150.00821744379970299419.825035769397502-8.68203858497106819.86407559226324-8.81313204267421421.046251732606027-8.31109348715200720.43269233852871-8.680831791344698
17513849.72977824536161619.8727132609065214.9149553653363,52.91916493746756377339221560167510381.015.96426909253809614.083802036553934-47.530699521812620.4708603535802704523.780048213690310.03256333246827125523.727838692607870.030124756018157582.2943933089856290.0046307528645993316.3219534915797390.00824851343420725220.45946817297725-7.16019546014345720.461854546852575-7.24221417903409622.998330331299535-6.73972978153553821.897871758596594-7.21260470769200444.779971644943460.0325484313070774144.4962773655858260.030555728261474645.0809373381893650.00714255510694007313.2379271949424170.01156215331553708319.772290464389595-7.84717326319009319.779190803289865-7.90882292141645922.135140402729373-7.13175080557387821.09545002936174-7.64790104118410550.610276038969620.01450593583285808650.335486423693690.0138515927890737025.3632260089096730.003696912408236274314.8949113045053850.00601978918805966919.639403234880845-8.85704353457783219.645314325424007-8.90123512321630822.076434753185204-7.904715493052254520.967405196392328-8.48408030951125
9482134.22578395710842975.9241914099766214.9160350651248,52.93737607078788420762203291730399111.014.4251549843533789.68810866411639263.57916700344990.74090295189763817.579878650712850.0269109830260276817.0574855066752630.021695766335087664.5930222517647960.0060609357030391698.1403942633359280.00865654775047655320.787460317681738-7.039376532484620.82021247264302-7.24023013297276822.244753626388103-7.20032897059360721.623386401092525-7.43440574186162135.796903022080010.02844637073576450334.8180981892626560.02458253310114512411.6953776261955460.01155676054558343718.790308133639510.01437602770261844620.015386362063573-7.75040892202431620.045487385779744-7.87871232174228321.2299643772781-7.51401767360486920.715165245112964-7.791567856333065551.25188470672450.01577029004693031349.91319981303220.014410070436411519.4437415953709730.00976384784927151730.789791860328920.0110145236086655519.625725398720498-8.78000943037523419.654461469213768-8.84918668128430520.67805039818143-8.24844216917776620.178983117317852-8.616490921519103
20091954.0848576295123388.545361033297214.91330171620754,52.94074700693218518812022333934429011.014.90922195423029210.10854632509445720.6579208256146420.735056096608601710.090487133700520.025555444881329.7556934434289740.0210087807326624323.60723690994937130.0054342689465528755.6253551892082780.00759903717796480121.390219667566463-6.493818007167150521.42685463782867-6.66947884847373622.507063334732944-7.05671826544326122.02462512666677-7.17494413218140821.1042412357444850.02550994418561458620.4565536976249650.0218964128896323827.694286648400860.0090351162578954112.1322517675167650.01141888950477702220.589075643724527-7.29546218628121420.622918824430155-7.42731027671213921.684579095534474-7.32686073349867821.190146464558907-7.56681527473796732.1809135822180540.01249407138675451331.196416366351870.01106159237085160310.3935973986182760.00522017805683630618.106373977053080.00670979389818495320.131004077248456-8.52765741478724620.16473822983424-8.62610255536460421.358085274522953-8.2482466070550220.75542128440663-8.578208038367503
17493889.98568762683361571.0494826433974214.9151090666597,52.9186458202309638253959151416288810.013.30745774299267311.02633443092811129.625732656513630.559866114512701119.509668390316620.0300899203866720219.202317990457970.025006477777386436.7082860774394950.00740830041122202110.5069415178335530.01002215611496325720.674375280829384-7.031245372294204520.691615856665383-7.214665826578877521.833471062539594-7.39343084704556221.34630921260184-7.55232302598089632.6062650025535350.0292906966060432.0142449146993040.0255281447907611611.173937197754530.0110638649666219918.00065342638540.01393732982684612820.116747365272843-7.617403278984502520.136641842380783-7.74667544212393521.27948443471843-7.51182290770977320.76177932414134-7.77861205702530641.98986580747420.0144003741443157241.246419959310610.0131477956686367915.70097613005630.00783307019759496625.5131714143677720.00926367364379397819.84213878322931-8.66229906513157719.86153435168754-8.74167931961506520.910183366633127-8.25552811549990120.383088882664875-8.600347153337044
.......................................................................................................................................................................................
18631033.58462759393572435.1291903129454214.93159110569908,52.941443516203041022104924272445282.06.09470765000732852.604613411647479324.395684946475490.90408315676070760.125476413147957380.0045182579196989540.147252853128352180.0094204987765787640.058044563764049330.0023748012722035420.12049003557383190.00479626239299968226.153594761202747-3.647386327320765625.97984070424018-3.05230309353709126.99059612102928-3.513872209827934726.197622168495297-3.54250160565664230.147145314755249970.0040874625556170940.17714226142317520.0084557474647795910.065838308311947690.00225908689775994450.139380531373330120.0043367755921382725.980633904219037-3.9204803515393325.779194538395053-3.35355327172857126.853803341672975-3.697993908231181226.039494710279552-3.80085530149379470.16984123587036770.00170954747591167690.242675106893440630.0035576579102228130.071638241910986470.00093868185260047650.164468264970689250.001802941575439186325.82489214658512-5.00377884847665525.437436926014815-4.600454139434086526.76213770200887-4.72070033626999225.859794722669758-4.912088486680949
1966721.22225340744492983.68670221368214.92861140426444,52.9463880447225171173229742992281.05.5657642972095053.371606007897061-41.4967647701128260.79563529983262490.11160096897987760.0044858786277472970.19015111821509810.0129877237876378750.042165922199666810.00238173614975969140.110198391639659970.00472579382029469426.28083008650714-3.532338747442450425.70225279076091-2.98564982147081327.33759599223483-3.179828079653862726.29456186006081-3.46484166064914940.106539448158504910.0040113953873515130.182537461674714630.0113697278360420360.044336733078342570.00214341279829961120.106424122462145940.004223173489675156526.331223892642214-3.6006664749622625.746619982415552-3.079609607246954527.283090776831408-3.340404157131588626.33239980545914-3.54575470531382830.117383047547148820.00170332915149629120.255966233772794540.0048555239369081720.044346525441187950.00089852470824821840.12079640998603860.001798654535931637426.22598654835853-4.61140875627747625.379543304038275-4.3252692947159627.282851003889096-4.25510262445816226.194864931390306-4.583813073208262
18132700.7722207463181958.4359030956539214.92167787788975,52.928286753555132693270819481974280.05.47516566902800952.7133103471874906-72.009727812952190.86856982657648470.119914876887652110.0046238433569669720.1552051358969640.0100104744343411260.071722027864496230.0024945592932914260.119188910790167440.00488959280199056926.202817335217446-3.57575321412237725.922734778976693-3.043991045588201826.76086859929756-3.683767870803752226.209410372535363-3.5110596348506660.174251041910502720.0041349390521645550.245556141717148180.0087326637582189560.100977410914683790.0022714192384853720.171984733671477620.00440281781017782425.7970615410429-4.08722902097951825.424622997645745-3.66045125815893926.389439422284166-4.14396959946754625.81127525448242-4.0068430425513380.15012188173982760.00166725972667336460.225928149722535140.00364376436664256950.07978227846659710.00092098604688213040.150458155123486450.001769110741591318325.958890001133305-4.89809361205879925.51507413598458-4.4984211769708626.64523391256543-4.85659513011477125.95646066920916-4.836843569272366
4622461.9428973677261399.0566144663642214.92970729322886,52.926782860239082454247113901409280.04.15514423820598563.8524813067274786-74.402110719027530.374666955344710230.17918792720555150.005010483320802450.210666307504597750.0090655910939056640.084916054837261710.0027039905014693150.18322311733050580.005290246610399315525.766728135744927-3.91351564381382825.59101229214823-3.46124245518759726.57752547792218-3.7764957796910125.742549330277267-3.87966556423802530.160515721054372980.0071357409469783310.186212113014930130.011955354550008110.083152502451970990.0044828536770624010.165975378431064220.0074483400685935325.88620606495768-3.42742067458269925.724980179433423-3.04867474824105526.600311690683355-3.22781198098962125.849890831174125-3.41762249251636070.156638591313128820.00171031721401959660.19374878119816950.00310196079704416460.075955142839186530.00093114765091898220.161247429925874760.00181158396844972625.912753078083792-4.916346018401225.681902551997023-4.50625187820515326.69860703800976-4.79207589039556225.881267996815545-4.8857158570738495
670974.94785924914072060.843736686114214.93605866920348,52.9397947446257996698420512071279.03.90511409180692943.186772503573259-54.302494097413690.57797981235380680.114527723312238180.00478540314361453060.132590420222687740.0070158916111602630.083995629982641110.00260586431490593540.116633874377125610.00506456831498880726.25272343125438-3.491924238662819326.09371963241331-3.247055457334690626.5893582706725-3.80393396787135626.232938244043122-3.4518564680372580.210676896398244550.00417841365560889240.234437924328451150.00607959998921670.143212448310666040.002374534969895530.21270481147982280.00440251591257581125.59095772026907-4.27783664846318325.474930331098108-3.9931792984555326.010048076619565-4.46885994309001425.58055671511738-4.2324338670313750.23060801850079360.00175578135531395670.25704398019116920.00255702138123523440.145420247082944650.0010046663085193730.23473102361001260.001857874643160739325.492813989675028-5.30424509242680625.374981406455706-5.01643001044723125.99343780445889-5.40898284345052425.473573768416895-5.26244488614722
2122660.0878462053954144.34451557572214.91693026935283,52.9530419918126865166841344157277.05.4543259973625132.918058782324085355.260191050584420.84485271126331040.102071474887921520.0044822152704000470.3821143910100170.013979951438507720.0552579755140685040.0024256435375079890.10114783793142310.00476553742300962626.377739024548646-3.440189459759723324.944516514529653-3.630732441118028427.044012575927276-3.440563802383956626.387608489956317-3.36709731112831050.138650667255372450.0040220455266535280.44732227244480170.0122669789635871630.076808264697901310.00219145544919071730.138298024270146260.00428105979084977926.045195090104787-3.874734842041690724.773448693764283-3.9340806065270426.686480116601615-3.892232266646387626.047960060476942-3.8062613176305750.120985220023089040.00172230962198227640.33046017112312540.0052359133182700480.060555667722989110.00094438140807871440.12183370312411720.00182365950112666726.193169203433552-4.631899962748148525.102202191457938-4.51738456893149926.94461300679262-4.534320085338972526.185581388122873-4.578190564911856
7211301.24402042397562193.4780971540504214.93173831037686,52.93843978842371292131021852202275.04.004356311506943.112305339884491-46.518701742393150.62921675712676060.1325621234030520.0079513294622302060.130454442294799970.0127203149790806930.106109519412150310.004463352560128170.130478941447079830.00844670699991002626.093951369821934-3.11819539288752726.111352818790195-2.62842147530380226.335614131134605-3.48496849202398726.111148938272393-3.04023750741907640.196644719105137170.0050712265074253080.219695103399969780.0077462862511188620.138433245181029980.0029815417987090760.197936992399490030.005384791345236881525.665794280195-3.999068054099707725.545449056493347-3.669439896557695826.046899000405467-4.1901348037335125.6586825825281-3.94253757898461470.209472691942868720.00177866243757307530.248904734728550430.00281957337798420.132360862049935250.00101026096844887020.213405513274313060.001877517837968873325.59718146529483-5.18676492103843525.40991710513455-4.87685438613498926.095601032491725-5.30157058440284425.576985912212304-5.148564319166033
1244025.1672750016087499.7676489967398214.92517747001168,52.9120045178760640184034489509274.04.0911250141223113.4446667037096526-78.3555832838150.53950062734358170.129861369418181350.0047570769675076010.16797469767726520.009527124789224940.081615651558781070.002605220199485510.13475676595822380.00509662229356160226.116300054258346-3.629410989833429525.836890329713942-3.175602406681268326.620566369849488-3.773938239063640726.07612355061181-3.59597635745669560.18161434233183510.0045204134657979010.23610627232994630.0086001187541618430.111633265068170190.00253064004729218040.184455884974540630.00478876673442359925.75212464408553-4.036623284680814525.467231188570643-3.63535215978169826.2805159318908-4.13574609869777425.735268710659252-3.9919998299905940.150084029259491520.0018751614261418580.206260643024554880.00365464241892520970.085021796548414370.00098965186239692730.154782869268439650.002010526648588393425.95916379859733-4.7717207782900625.613959081641088-4.397997932831188526.57617430609045-4.84768453712908525.925692767302625-4.7300448313878665
16462922.4683172105324891.6568203762286214.93092533581597,52.9211200546841829122934884902272.04.3155491984752373.027292497526654-28.7091098998549620.71268433946256660.150247613730477350.0047318567521870140.204388669803026280.0092755749802644160.086160344887902490.0025937095779133550.147360473886383240.00505907173644474325.95798104210533-3.78810645994783125.623857956989212-3.40597752913908826.561731428215147-3.835667139649283625.979047476318282-3.69742456243442770.30274587334412320.0042685181833803650.40665462258935820.0080577296919110970.214295965639256620.00252349039182813960.310130074372062970.00455477358745599125.197304419769377-4.64220399506509124.876935716385482-4.27883569207287825.572465012398943-4.83524148751053525.17114029175494-4.5985225965450490.177998921665170060.00174600386526435610.2583400169686090.0033232626627179430.114222759485027040.00098643107159905750.184822509485832730.001850689673699664225.773956571703675-5.03152859409897625.369520790581404-4.74044528053597326.255618380025304-5.16855095424085525.733112843204136-5.00937092159065

Let’s save our source catalog table to an Astropy ECSV file. The ECSV format is the recommended way to store Table data in a human-readable ASCII file.

We’ll use this file later in the Imviz tutorial.

tmp_tbl.write('photutils_cat.ecsv', overwrite=True)

Photometry types#

Segment fluxes (aka “Isophotal Fluxes”)#

The segment fluxes are simply the sources fluxes measured using the pixels within the source segments.

This is sometimes referred to as “isophotal flux” because source segments were defined as the pixels above a threshold flux level (i.e. isophote).

label = segm_deblended.data[1130, 3780]
idx = segm_deblended.get_index(label)
segm_deblended.polygons[idx]
../../../_images/e0a0b48add146fc4c586ff5070c07be671d36212175ddcab9a83933373b67975.svg
fig, ax = plt.subplots(ncols=4, figsize=(12, 3))
norm3 = simple_norm(f200w_cat.data[idx].value, 'sqrt', percent=98)
ax[0].imshow(f200w_cat.data[idx].value, norm=norm3)
ax[0].set_title('F200W Data')
norm3 = simple_norm(f200w_cat.error[idx].value, 'sqrt', percent=98)
ax[1].imshow(f200w_cat.error[idx].value, norm=norm3)
ax[1].set_title('F200W Error')
ax[2].imshow(f200w_cat.segment_ma[idx])
ax[2].set_title('Segment (detection image)')
norm3 = simple_norm(det_cat.data[idx].value, 'sqrt', percent=90)
ax[3].imshow(det_cat.data_ma[idx], norm=norm3)
_ = ax[3].set_title('Detection Image)')
fig, ax = plt.subplots(ncols=4, figsize=(12, 3))
norm3 = simple_norm(f200w_cat.data[idx].value, 'sqrt', percent=98)
ax[0].imshow(f200w_cat.data_ma[idx], norm=norm3)
ax[0].set_title('F200W Data')
norm3 = simple_norm(f200w_cat.error[idx].value, 'sqrt', percent=98)
ax[1].imshow(f200w_cat.error_ma[idx], norm=norm3)
ax[1].set_title('F200W Error')
ax[2].imshow(f200w_cat.segment_ma[idx])
ax[2].set_title('Segment (detection image)')
norm3 = simple_norm(det_cat.data[idx].value, 'sqrt', percent=90)
ax[3].imshow(det_cat.data_ma[idx], norm=norm3)
_ = ax[3].set_title('Detection Image)')

Aperture Photometry#

Aperture fluxes are measured in a circular aperture. We can use the plot_circular_apertures method to plot them. Recall we computed aperture fluxes above with r=5 and r=10 pixels.

fig, ax = plt.subplots(figsize=(8, 8))
norm = simple_norm(f200w['sci_bksub'].value, 'sqrt', percent=99)
ax.imshow(f200w['sci_bksub'].value, norm=norm)
_ = det_cat.plot_circular_apertures(radius=10, ax=ax, color='white')

Kron Photometry#

Kron photometry is measured in elliptical apertures. The shape of the elliptical apertures are different for each source and defined from the object’s shape.

The object’s shape is based on calculating image moments from the object’s profile.

The size of the aperture calculated from the Kron radius (measured from a first moment of the object’s profile) times a scale factor (which can be specified). It is defined as the radius within which a fixed fraction (around 90%) of the total light emitted by the object is enclosed.

Let’s use the plot_kron_apertures method to plot the Kron apertures.

fig, ax = plt.subplots(figsize=(8, 8))
ax.imshow(f200w['sci_bksub'].value, norm=norm)
_ = det_cat.plot_kron_apertures(ax=ax, color='white')

Let’s export the Kron apertures for the largest 500 sources. We’ll save the apertures to a DS9 region file.

idx = segm_deblended.get_index(tmp_tbl['label'])
objs = det_cat[idx]
objs
<photutils.segmentation.catalog.SourceCatalog>
Length: 500
labels: [1922 1757  320 1791 2140 ...  670 2122  721  124 1646]
fig, ax = plt.subplots(figsize=(8, 8))
ax.imshow(f200w['sci_bksub'].value, norm=norm)
_ = objs.plot_kron_apertures(ax=ax, color='white')

Here, we convert the Kron apertures to Astropy Region objects, which can be serialized to DS9 format.

We’ll use this file later in the Imviz tutorial.

regions = Regions([aperture_to_region(aper.to_sky(wcs=f356w['wcs'])) for aper in objs.kron_aperture])
regions[0:3]  # show the first 3 regions
<Regions([<EllipseSkyRegion(center=<SkyCoord (ICRS): (ra, dec) in deg
    (214.92110659, 52.93793595)>, width=4.684728206075569 arcsec, height=3.4754142265185948 arcsec, angle=-1.3715084074762762 rad)>, <EllipseSkyRegion(center=<SkyCoord (ICRS): (ra, dec) in deg
    (214.9152941, 52.92080683)>, width=5.947715306825106 arcsec, height=1.6448272422015742 arcsec, angle=0.20554718357823298 rad)>, <EllipseSkyRegion(center=<SkyCoord (ICRS): (ra, dec) in deg
    (214.92093216, 52.91691942)>, width=3.122354338355514 arcsec, height=2.3661185452014624 arcsec, angle=0.16841946780471906 rad)>])>
regions.write('kron_apertures.reg', overwrite=True)

Data Exploration#

As an example, let’s look for Lyman-break F115W dropout galaxies in our multiband catalog.

We’ll use a color-color selection criteria, where the first color is F115W - F200W (the dropout color) and the second is F200W - F356W.

f115w_f200w = mult_tbl['f115w_kron_mag'] - mult_tbl['f200w_kron_mag']
f200w_f356w = mult_tbl['f200w_kron_mag'] - mult_tbl['f356w_kron_mag']

f115w_segm_snr = mult_tbl['f115w_segment_flux'] / mult_tbl['f115w_segment_fluxerr']
f115w_snr = mult_tbl['f115w_kron_flux'] / mult_tbl['f115w_kron_fluxerr']
f200w_snr = mult_tbl['f200w_kron_flux'] / mult_tbl['f200w_kron_fluxerr']
f356w_snr = mult_tbl['f356w_kron_flux'] / mult_tbl['f356w_kron_fluxerr']

col_t1 = 1.5
col_t2 = 0.0
mask = ((f115w_f200w.value > col_t1) & (f200w_f356w.value < col_t2)
        & (f115w_segm_snr < 2) & (f115w_snr < 2) & (f200w_snr > 5) & (f356w_snr > 5))


fig, ax = plt.subplots()
ax.scatter(f200w_f356w, f115w_f200w, s=1)
ax.scatter(f200w_f356w[mask], f115w_f200w[mask], s=4, color='red')
ax.set_xlabel('F200W – F356W')
ax.set_ylabel('F115W – F200W')
ax.plot((-5, col_t2), (col_t1, col_t1), color='red', ls='dashed')
ax.plot((col_t2, col_t2), (col_t1, 10), color='red', ls='dashed')
ax.set_xlim(-5, 7)
ax.set_ylim(-6, 8)
_ = ax.set_title('Color-color selection')
hiz_tbl = mult_tbl[mask]
hiz_tbl
QTable length=5
labelxcentroidycentroidsky_centroidbbox_xminbbox_xmaxbbox_yminbbox_ymaxareasemimajor_sigmasemiminor_sigmaorientationeccentricityf115w_segment_fluxf115w_segment_fluxerrf115w_kron_fluxf115w_kron_fluxerrf115w_circ_rad5_fluxf115w_circ_rad5_fluxerrf115w_circ_rad10_fluxf115w_circ_rad10_fluxerrf115w_segment_magf115w_segment_magerrf115w_kron_magf115w_kron_magerrf115w_circ_rad5_magf115w_circ_rad5_magerrf115w_circ_rad10_magf115w_circ_rad10_magerrf200w_segment_fluxf200w_segment_fluxerrf200w_kron_fluxf200w_kron_fluxerrf200w_circ_rad5_fluxf200w_circ_rad5_fluxerrf200w_circ_rad10_fluxf200w_circ_rad10_fluxerrf200w_segment_magf200w_segment_magerrf200w_kron_magf200w_kron_magerrf200w_circ_rad5_magf200w_circ_rad5_magerrf200w_circ_rad10_magf200w_circ_rad10_magerrf356w_segment_fluxf356w_segment_fluxerrf356w_kron_fluxf356w_kron_fluxerrf356w_circ_rad5_fluxf356w_circ_rad5_fluxerrf356w_circ_rad10_fluxf356w_circ_rad10_fluxerrf356w_segment_magf356w_segment_magerrf356w_kron_magf356w_kron_magerrf356w_circ_rad5_magf356w_circ_rad5_magerrf356w_circ_rad10_magf356w_circ_rad10_magerr
deg,degpix2pixpixdeguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmaguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmaguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmag
int32float64float64SkyCoordint64int64int64int64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64
616401.19550590719951876.4877753657443214.94314150986483,52.9424407821729943994041874187928.01.44145315510549591.409034122127981279.23981215692240.21089148672542876-0.00095152588398756650.001428473275154829-0.0065199256812457990.005077226330820457-0.00321235506191488950.0024014290033907127-0.0060026645297615690.00478829411802379430.260244369506836nan28.883358701680752nan29.696250631977886nan28.94697296377662nan0.0140174303684064780.00129497447051107880.043417640735081240.004511952315104780.0232890368377580880.00216874715075785670.041610188278918530.004251772514169767528.5333289814555-2.68195549008598527.30583444826418-2.56559810117667427.98212118041259-2.674029188922715727.352000797737627-2.58220649658969760.00767754711714257850.00050593452760949730.0270075125344489760.00179864600641908670.0143099851404636810.00084191127293364490.0255723783554595440.001693423671576686729.186943774023618-3.022109484264982327.821288534191698-3.01134912207767628.51090204303207-3.138001669919714427.880572196494786-3.017131626902192
18222966.13448516637252034.2413113622565214.91850323172991,52.92701134774817296329702030203960.02.48951051920270141.8060176375926853-66.174490842858790.68827392585057970.0034158989167557640.0020931400358676910.0035225558279546280.0080773611244366550.0040014439060259480.0023984409418828960.0071390766477886610.00477285119266031430.066237473576376-1.05069389632997330.032855287895345-0.3929634503322998529.89445816713882-1.065607826999590129.26589488817223-0.99301039925556960.016478060786049450.0019270342309027910.067836286243794370.0070915188626633680.0198679684666919660.00217688207006882070.044393213280670290.00419293791007389628.35773474850751-2.45012158632059926.82134483981966-2.559759370986012328.154616344909336-2.513680640919897527.281708546562644-2.65998516487746750.0147860692511502770.00080581504153087740.06186732924468570.00308533478028395370.0178378102169849950.00092562029648201860.038965757283537910.001832358435897911628.475368160245143-3.216658083634799326.92134657916492-3.30823668351810828.271646152513384-3.267198435610112527.423292197518375-3.3690741837736757
20501112.59745689392573613.731159057536214.91846711453576,52.947304904449375111011163612361514.01.6637968228104110.7171950573725432-18.9113207883365570.90232361087916390.00110290992247945050.00096071138978004460.0052596620174963430.004543657959517810.00192774978971949890.00230859456273833660.0168792002515879520.00468149053716292531.293649890176944-0.830092662261619829.597605406085393-0.834918868164816830.68737333894067-0.659108965396594428.331620336004345-1.65819630336365530.0029761775832979190.00092152715660631660.049175866260453470.0043579817477491470.0163757766501282480.0022303093465082170.053411920474156940.00449885483993938430.215852896906956-1.56575199799140727.170619952501166-2.72335763641781928.364495234162348-2.303224798853040327.08090451574303-2.77414349724177360.0031754741542403340.00036552429082803430.0284039587337129860.00168394548645310160.0101157071570862730.00085365660071076440.0293902058676936940.001718225256383744830.145478544280053-2.465523732156109627.76655281745347-3.13015977431216928.887509379196494-2.77224558926473527.72949342960494-3.1444950265393286
2094279.427217429013353906.0357183500323214.92285230495403,52.9541721212956842752843904390843.02.59135280026410661.2470181273581454.8076756253059850.87659830489896720.00119752718509658330.0017749716062098742-0.0029286628927670760.00712363431228187850.0010053583806207030.0023968218416314394-0.00232579360249307930.00480709500853115731.204286547177176-0.559825691380920628.515670968082404nan31.394197743664265-0.380303904825917428.94271824656635nan0.0122971592187169270.0015192308928817510.07560572706797170.0061784948783780160.0163175113476037480.0020360595725433370.040312063702855920.00416986286582597928.675488009638887-2.396927010542974326.70361326456007-2.804467088428623428.368365191531705-2.387325233877792427.38641242090234-2.57015454062350820.0092060724635467920.00062942359363660220.055199497233397470.0025185508276166970.0134270117113871620.00084230621256607990.02970839429984990.00168125564071100928.989814027441017-2.98463312606153327.045162194705785-3.400402288914668328.58005158070009-3.07233302585288727.717802051537607-3.1778817950353218
2101922.84820648370243908.8543242193277214.91705750417216,52.95010264310915591193639053912167.06.9821498383174791.76113480270446780.33366812486586220.96766629866495530.0005123345553004030.0033893389627337456-0.0294140412062027930.010624864016645653-2.4691226962951828e-050.0023222670564890567-0.0035182976441444020.00468309445069128532.12611587846169-0.152839804463266528.08161655984555nan29.732644607464486nan28.97109271856273nan0.0529349967181896860.00307786767370998860.215034078573094920.0094013326709368110.0212985858736242580.0020705962764032620.062954504752051040.00415408717784970327.090642773773244-3.15009459473376525.568731769240955-3.444754899872242628.079123076810372-2.631370212508288326.90243297052409-3.020756300358290.050768644496184320.00125552667304873470.153596897765700570.00390493354874997140.022737530420861060.00086719139508835640.055731870351664360.001706954423316535827.136011079851997-4.04344806597940725.934043889420828-4.01417990008194428.00814176719847-3.587209816916790727.034740954633996-3.817459050265213
shape = (51, 51)
labels = hiz_tbl['label']
fig, ax = plt.subplots(ncols=5, nrows=len(hiz_tbl), figsize=(8, 8))
for i, label in enumerate(labels):
    idx = segm_deblended.get_index(label)
    cimg1 = f115w_cat[idx].make_cutouts(shape)
    cimg2 = f200w_cat[idx].make_cutouts(shape)
    cimg3 = f356w_cat[idx].make_cutouts(shape)
    cimg4 = det_cat[idx].make_cutouts(shape)
    cimg5 = CutoutImage(segm_deblended.data, (det_cat[idx].ycentroid, det_cat[idx].xcentroid), shape)
    
    snorm = SimpleNorm('sqrt', percent=99)
    norm = snorm(cimg4.data)
    ax[i, 0].imshow(cimg1.data, norm=norm)
    idstr = f'ID = {label}'
    ax[i, 0].set_ylabel(idstr)
    ax[i, 0].set_title('F115W')
    ax[i, 1].imshow(cimg2.data, norm=norm)
    ax[i, 1].set_title('F200W')
    ax[i, 2].imshow(cimg3.data, norm=norm)
    ax[i, 2].set_title('F356W')
    ax[i, 3].imshow(cimg4.data, norm=norm)
    ax[i, 3].set_title('Detection')
    det_cat[idx].plot_kron_apertures(ax=ax[i, 3], origin=cimg4.xyorigin, color='white')
    ax[i, 4].imshow(cimg5, cmap=segm_deblended.cmap)
    ax[i, 4].set_title('Segment')

The best candidate shown above was reported by the CEERS team (Finkelstein et al. 2022) as a high-redshift galaxy candidate called “Maisie’s Galaxy” with a photometric redshift of z = 11.8.

The CEERS team later obtain a NIRSpec spectrum of Maisie’s Galaxy (see Arrabal Haro et al. 2023) and measured a spectroscopic redshift of z = 11.416.

# label = 616  # on linux
label = segm_deblended.data[1876, 401]
idx = np.where(hiz_tbl['label'] == label)[0][0]
obj = hiz_tbl[idx]
obj
Row index=0
labelxcentroidycentroidsky_centroidbbox_xminbbox_xmaxbbox_yminbbox_ymaxareasemimajor_sigmasemiminor_sigmaorientationeccentricityf115w_segment_fluxf115w_segment_fluxerrf115w_kron_fluxf115w_kron_fluxerrf115w_circ_rad5_fluxf115w_circ_rad5_fluxerrf115w_circ_rad10_fluxf115w_circ_rad10_fluxerrf115w_segment_magf115w_segment_magerrf115w_kron_magf115w_kron_magerrf115w_circ_rad5_magf115w_circ_rad5_magerrf115w_circ_rad10_magf115w_circ_rad10_magerrf200w_segment_fluxf200w_segment_fluxerrf200w_kron_fluxf200w_kron_fluxerrf200w_circ_rad5_fluxf200w_circ_rad5_fluxerrf200w_circ_rad10_fluxf200w_circ_rad10_fluxerrf200w_segment_magf200w_segment_magerrf200w_kron_magf200w_kron_magerrf200w_circ_rad5_magf200w_circ_rad5_magerrf200w_circ_rad10_magf200w_circ_rad10_magerrf356w_segment_fluxf356w_segment_fluxerrf356w_kron_fluxf356w_kron_fluxerrf356w_circ_rad5_fluxf356w_circ_rad5_fluxerrf356w_circ_rad10_fluxf356w_circ_rad10_fluxerrf356w_segment_magf356w_segment_magerrf356w_kron_magf356w_kron_magerrf356w_circ_rad5_magf356w_circ_rad5_magerrf356w_circ_rad10_magf356w_circ_rad10_magerr
deg,degpix2pixpixdeguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmaguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmaguJyuJyuJyuJyuJyuJyuJyuJymagmagmagmagmagmagmagmag
int32float64float64SkyCoordint64int64int64int64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64float64float32float64float64float64float64float64float64float64float64float64float64float64float64float64float64
616401.19550590719951876.4877753657443214.94314150986483,52.9424407821729943994041874187928.01.44145315510549591.409034122127981279.23981215692240.21089148672542876-0.00095152588398756650.001428473275154829-0.0065199256812457990.005077226330820457-0.00321235506191488950.0024014290033907127-0.0060026645297615690.00478829411802379430.260244369506836nan28.883358701680752nan29.696250631977886nan28.94697296377662nan0.0140174303684064780.00129497447051107880.043417640735081240.004511952315104780.0232890368377580880.00216874715075785670.041610188278918530.004251772514169767528.5333289814555-2.68195549008598527.30583444826418-2.56559810117667427.98212118041259-2.674029188922715727.352000797737627-2.58220649658969760.00767754711714257850.00050593452760949730.0270075125344489760.00179864600641908670.0143099851404636810.00084191127293364490.0255723783554595440.001693423671576686729.186943774023618-3.022109484264982327.821288534191698-3.01134912207767628.51090204303207-3.138001669919714427.880572196494786-3.017131626902192
fluxes = u.Quantity((obj['f115w_kron_flux'], obj['f200w_kron_flux'], obj['f356w_kron_flux']))
fluxes.to(u.nJy)
\[[-6.5199257,~43.417641,~27.007513] \; \mathrm{nJy}\]
# measured fluxes from Finkelstein et al. 2022
ceers_fluxes = np.array((-7.52, 45.11, 26.29)) * u.nJy
ceers_fluxes
\[[-7.52,~45.11,~26.29] \; \mathrm{nJy}\]
fluxes.to(u.ABmag)
/opt/hostedtoolcache/Python/3.11.13/x64/lib/python3.11/site-packages/astropy/units/function/logarithmic.py:67: RuntimeWarning: invalid value encountered in log10
  return dex.to(self._function_unit, np.log10(x))
\[[{\rm NaN},~27.305834,~27.821289] \; \mathrm{mag}$$\mathrm{\left( \mathrm{AB} \right)}\]
ceers_fluxes.to(u.ABmag)
\[[{\rm NaN},~27.264318,~27.850524] \; \mathrm{mag}$$\mathrm{\left( \mathrm{AB} \right)}\]
Space Telescope Logo