Petrosian

class petrofit.petrosian.Petrosian(r_list, area_list, flux_list, area_err=None, flux_err=None, eta=0.2, epsilon=2.0, epsilon_fraction=0.99, total_flux_fraction=0.99, verbose=False, interp_kind='cubic', interp_num=5000)[source]

Bases: object

Class that computes and plots Petrosian properties.

Petrosian properties.

Parameters:
r_listnumpy.array

Array of radii in pixels.

area_listnumpy.array

Array of aperture areas.

flux_listnumpy.array

Array of photometric flux values.

area_errnumpy.array

Array of errors in the aperture areas. If None and flux_err is provided, errors are calculated with area_err set to zero.

flux_errnumpy.array

Array of errors in the flux values. If None and area_err is provided, errors are calculated with flux_err is set to zero.

etafloat, default=0.2

Eta is the petrosian value which defines the r_petrosian.

epsilonfloat

Epsilon value (used to determine r_total_flux). N.B: r_total_flux = r_petrosian * epsilon

epsilon_fraction: float, default=0.99

Fraction of total flux recovered by r_petrosian * epsilon.

total_flux_fractionfloat, default=0.99

Fraction of Sersic flux that defines the Petrosian total flux. Sersic total flux is the flux at infinity, thus a smaller fraction must be used to define the total flux when analysing images. total_flux_fraction can also be adjusted if the image has low signal-to-noise or if the profile extends too far out (for example profiles with high Sersic indices).

verbosebool

Prints info

interp_kindstr or int, optional

Specifies the kind of interpolation used on the curve of growth (i.e r_list vs flux_list)

interp_numint

Number of interpolation function sampling radii.

Attributes Summary

area_err

area_list

c2080

c2080 = 5 * np.log10(r_80 / r_20)

c5090

c5090 = 5 * np.log10(r_90 / r_50)

epsilon

Epsilon value (used to determine r_total_flux).

epsilon_fraction

Fraction of total flux recovered by epsilon

eta

Eta is the Petrosian value which defines the r_petrosian

flux_err

flux_list

half_flux

Returns the flux at r_half_flux or np.nan if out of range

half_flux_err

Returns the photometric uncertainty at r_half_flux or np.nan if out of range.

r_50

r_50_err

r_epsilon

r_epsilon = r_petrosian * epsilon

r_epsilon_in_range

True if r_epsilon is in input radii list

r_half_light

R_e or Radius containing half of the total Petrosian flux

r_half_light_err

Error estimate on r_e

r_list

r_petrosian

The Petrosian radius is defined as the radius at which the Petrosian profile equals eta.

r_petrosian_err

Estimated 1-sigma r_petrosian Error.

r_total_flux

The radius containing the total flux

r_total_flux_err

The radius containing the total flux

total_flux

Returns the flux at r_total_flux or np.nan if out of range

total_flux_err

Returns the photometric uncertainty at r_total_flux or np.nan if out of range.

total_flux_fraction

Fraction of Sersic flux that defines the Petrosian total flux.

Methods Summary

concentration_index([fraction_1, fraction_2])

Calculates Petrosian concentration index.

fraction_flux_to_r([fraction])

Given a fraction, compute the radius containing that fraction of the Petrosian total flux

fraction_flux_to_r_err([fraction])

Given a fraction, compute the radius containing that fraction of the Petrosian total flux

imshow([position, elong, theta, color, lw])

Make 2D plots of elliptical apertures with radii of r_half_light, r_total_flux, r_20 and r_80.

plot([plot_r, title, radius_unit, ax, ...])

Plots the Petrosian profile.

plot_cog([plot_r, title, radius_unit, ...])

Plots the Curve of Growth (COG) for the Petrosian profile.

r_half_light_arcsec(wcs)

Given a wcs, compute the radius containing half of the total Petrosian flux in arcsec.

r_total_flux_arcsec(wcs)

Given a wcs, compute the radius containing the total Petrosian flux in arcsec.

Attributes Documentation

area_err
area_list
c2080

c2080 = 5 * np.log10(r_80 / r_20)

c5090

c5090 = 5 * np.log10(r_90 / r_50)

epsilon

Epsilon value (used to determine r_total_flux). N.B: r_total_flux = r_petrosian * epsilon

epsilon_fraction

Fraction of total flux recovered by epsilon

eta

Eta is the Petrosian value which defines the r_petrosian

flux_err
flux_list
half_flux

Returns the flux at r_half_flux or np.nan if out of range

half_flux_err

Returns the photometric uncertainty at r_half_flux or np.nan if out of range. Note that this includes errors in estimating r_half_flux.

r_50
r_50_err
r_epsilon

r_epsilon = r_petrosian * epsilon

r_epsilon_in_range

True if r_epsilon is in input radii list

r_half_light

R_e or Radius containing half of the total Petrosian flux

r_half_light_err

Error estimate on r_e

r_list
r_petrosian

The Petrosian radius is defined as the radius at which the Petrosian profile equals eta.

r_petrosian_err

Estimated 1-sigma r_petrosian Error.

r_total_flux

The radius containing the total flux

r_total_flux_err

The radius containing the total flux

total_flux

Returns the flux at r_total_flux or np.nan if out of range

total_flux_err

Returns the photometric uncertainty at r_total_flux or np.nan if out of range. Note that this includes errors in estimating r_total_flux.

total_flux_fraction

Fraction of Sersic flux that defines the Petrosian total flux. Sersic total flux is the flux at infinity, thus a smaller fraction must be used to define the total flux when analysing images. total_flux_fraction can also be adjusted if the image has low signal-to-noise or if the profile extends too far out (for example profiles with high Sersic indices).

Methods Documentation

concentration_index(fraction_1=0.2, fraction_2=0.8)[source]

Calculates Petrosian concentration index.

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

Parameters:
fraction_1float

Fraction of total light enclosed by the radius in the denominator.

fraction_2float

Fraction of total light enclosed by the radius in the numerator.

Returns:
r_fraction_1, r_fraction_2, concentration_index
  • r_fraction_1float or np.nan

    Radius containing fraction_1 of the total flux.

  • r_fraction_2: float or np.nan

    Radius containing fraction_2 of the total flux.

  • concentration_indexfloat or np.nan

    Concentration index

fraction_flux_to_r(fraction=0.5)[source]

Given a fraction, compute the radius containing that fraction of the Petrosian total flux

fraction_flux_to_r_err(fraction=0.5)[source]

Given a fraction, compute the radius containing that fraction of the Petrosian total flux

imshow(position=(0, 0), elong=1.0, theta=0.0, color=None, lw=None)[source]

Make 2D plots of elliptical apertures with radii of r_half_light, r_total_flux, r_20 and r_80.

Parameters:
positiontuple

(x, y) center of the apertures.

elongfloat

Elongation of the aperture.

thetafloat

The orientation of the aperture in rad.

colorstr

Color override that will change the color of the apertures in the plot.

lwfloat

Line width (thickness) of the plotted apertures.

plot(plot_r=True, title='Petrosian Profile', radius_unit='pix', ax=None, color='tab:blue', err_alpha=0.2, err_capsize=3, show_legend=True, legend_fontsize=None, ax_fontsize=None, tick_fontsize=None)[source]

Plots the Petrosian profile.

Parameters:
plot_rbool, optional

If True, plots the total flux and half-light radii. Default is True.

titlestr, optional

Title for the plot. Default is ‘Petrosian Profile’.

radius_unitstr, optional

Unit for the radius. Default is ‘pix’.

axmatplotlib.axis, optional

Matplotlib axis object to plot on. If None, creates a new axis.

colorstring

Matplotlib color for profile.

err_alphafloat, optional

Transparency for the error region. Default is 0.2.

err_capsizeint, optional

Cap size for the error bars. Default is 3.

show_legendbool, optional

If True, displays the legend. Default is True.

legend_fontsizeint or float, optional

Font size for the legend. If None, uses default font size.

ax_fontsizeint or float, optional

Font size for the axis labels. If None, uses default font size.

tick_fontsizeint or float, optional

Font size for the tick labels. If None, uses default font size.

Returns:
axmatplotlib.axis

Matplotlib axis object with the plot.

plot_cog(plot_r=True, title='Curve of Growth', radius_unit='pix', flux_unit='', ax=None, color='tab:blue', err_alpha=0.2, err_capsize=3, show_legend=True, legend_fontsize=None, ax_fontsize=None, tick_fontsize=None)[source]

Plots the Curve of Growth (COG) for the Petrosian profile.

Parameters:
plot_rbool, optional

If True, plots radii of interest including Petrosian radius. Default is True.

titlestr, optional

Title for the plot. Default is ‘Curve of Growth’.

radius_unitstr, optional

Unit for the radius. Default is ‘pix’.

flux_unitstr, optional

Unit for the cumulative flux. Default is ‘’

axmatplotlib.axis, optional

Matplotlib axis object to plot on. If None, creates a new axis.

colorstring

Matplotlib color for profile.

err_alphafloat, optional

Transparency for the error region. Default is 0.2.

err_capsizeint, optional

Cap size for the error bars. Default is 3.

show_legendbool, optional

If True, displays the legend. Default is True.

legend_fontsizeint or float, optional

Font size for the legend. If None, uses default font size.

ax_fontsizeint or float, optional

Font size for the axis labels. If None, uses default font size.

tick_fontsizeint or float, optional

Font size for the tick labels. If None, uses default font size.

Returns:
axmatplotlib.axis

Matplotlib axis object with the plot.

r_half_light_arcsec(wcs)[source]

Given a wcs, compute the radius containing half of the total Petrosian flux in arcsec.

r_total_flux_arcsec(wcs)[source]

Given a wcs, compute the radius containing the total Petrosian flux in arcsec.