LSST Applications g0f08755f38+9c285cab97,g1635faa6d4+13f3999e92,g1653933729+a8ce1bb630,g1a0ca8cf93+bf6eb00ceb,g28da252d5a+0829b12dee,g29321ee8c0+5700dc9eac,g2bbee38e9b+9634bc57db,g2bc492864f+9634bc57db,g2cdde0e794+c2c89b37c4,g3156d2b45e+41e33cbcdc,g347aa1857d+9634bc57db,g35bb328faa+a8ce1bb630,g3a166c0a6a+9634bc57db,g3e281a1b8c+9f2c4e2fc3,g414038480c+077ccc18e7,g41af890bb2+fde0dd39b6,g5fbc88fb19+17cd334064,g781aacb6e4+a8ce1bb630,g80478fca09+55a9465950,g82479be7b0+d730eedb7d,g858d7b2824+9c285cab97,g9125e01d80+a8ce1bb630,g9726552aa6+10f999ec6a,ga5288a1d22+2a84bb7594,gacf8899fa4+c69c5206e8,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gc28159a63d+9634bc57db,gcf0d15dbbd+4b7d09cae4,gda3e153d99+9c285cab97,gda6a2b7d83+4b7d09cae4,gdaeeff99f8+1711a396fd,ge2409df99d+5e831397f4,ge79ae78c31+9634bc57db,gf0baf85859+147a0692ba,gf3967379c6+41c94011de,gf3fb38a9a8+8f07a9901b,gfb92a5be7c+9c285cab97,w.2024.46
LSST Data Management Base Package
Loading...
Searching...
No Matches
Classes | Functions
lsst.pipe.tasks.computeExposureSummaryStats Namespace Reference

Classes

class  ComputeExposureSummaryStatsConfig
 
class  ComputeExposureSummaryStatsTask
 

Functions

 maximum_nearest_psf_distance (image_mask, psf_cat, sampling=8, bad_mask_bits=["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"])
 
 compute_psf_image_deltas (image_mask, image_psf, sampling=96, ap_radius_pix=3.0, bad_mask_bits=["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"])
 
 compute_ap_corr_sigma_scaled_delta (image_mask, image_ap_corr_field, psfSigma, sampling=96, bad_mask_bits=["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"])
 
 compute_magnitude_limit (psfArea, skyBg, zeroPoint, readNoise, gain, snr)
 

Function Documentation

◆ compute_ap_corr_sigma_scaled_delta()

lsst.pipe.tasks.computeExposureSummaryStats.compute_ap_corr_sigma_scaled_delta ( image_mask,
image_ap_corr_field,
psfSigma,
sampling = 96,
bad_mask_bits = ["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"] )
Compute the delta between the maximum and minimum aperture correction
values scaled (divided) by ``psfSigma`` for the given field representation,
``image_ap_corr_field`` evaluated on a grid of points lying in the
unmasked region of the image.

Parameters
----------
image_mask : `lsst.afw.image.Mask`
    The mask plane associated with the exposure.
image_ap_corr_field : `lsst.afw.math.ChebyshevBoundedField`
    The ChebyshevBoundedField representation of the aperture correction
    of interest for the exposure.
psfSigma : `float`
    The PSF model second-moments determinant radius (center of chip)
    in pixels.
sampling : `int`, optional
    Sampling rate in each dimension to create the grid of points at which
    to evaluate ``image_psf``s trace radius value. The tradeoff is between
    adequate sampling versus speed.
bad_mask_bits : `list` [`str`], optional
    Mask bits required to be absent for a pixel to be considered
    "unmasked".

Returns
-------
ap_corr_sigma_scaled_delta : `float`
    The delta between the maximum and minimum of the (multiplicative)
    aperture correction values scaled (divided) by ``psfSigma`` evaluated
    on the x,y-grid subsampled on the unmasked detector pixels by a factor
    of ``sampling``.  If the aperture correction evaluates to NaN on any
    of the grid points, this is set to NaN.

Definition at line 845 of file computeExposureSummaryStats.py.

851):
852 """Compute the delta between the maximum and minimum aperture correction
853 values scaled (divided) by ``psfSigma`` for the given field representation,
854 ``image_ap_corr_field`` evaluated on a grid of points lying in the
855 unmasked region of the image.
856
857 Parameters
858 ----------
859 image_mask : `lsst.afw.image.Mask`
860 The mask plane associated with the exposure.
861 image_ap_corr_field : `lsst.afw.math.ChebyshevBoundedField`
862 The ChebyshevBoundedField representation of the aperture correction
863 of interest for the exposure.
864 psfSigma : `float`
865 The PSF model second-moments determinant radius (center of chip)
866 in pixels.
867 sampling : `int`, optional
868 Sampling rate in each dimension to create the grid of points at which
869 to evaluate ``image_psf``s trace radius value. The tradeoff is between
870 adequate sampling versus speed.
871 bad_mask_bits : `list` [`str`], optional
872 Mask bits required to be absent for a pixel to be considered
873 "unmasked".
874
875 Returns
876 -------
877 ap_corr_sigma_scaled_delta : `float`
878 The delta between the maximum and minimum of the (multiplicative)
879 aperture correction values scaled (divided) by ``psfSigma`` evaluated
880 on the x,y-grid subsampled on the unmasked detector pixels by a factor
881 of ``sampling``. If the aperture correction evaluates to NaN on any
882 of the grid points, this is set to NaN.
883 """
884 mask_arr = image_mask.array[::sampling, ::sampling]
885 bitmask = image_mask.getPlaneBitMask(bad_mask_bits)
886 good = ((mask_arr & bitmask) == 0)
887
888 x = np.arange(good.shape[1], dtype=np.float64) * sampling
889 y = np.arange(good.shape[0], dtype=np.float64) * sampling
890 xx, yy = np.meshgrid(x, y)
891
892 ap_corr = image_ap_corr_field.evaluate(xx.ravel(), yy.ravel()).reshape(xx.shape)
893 ap_corr_good = ap_corr[good]
894 if ~np.isfinite(ap_corr_good).all():
895 ap_corr_sigma_scaled_delta = float("nan")
896 else:
897 ap_corr_sigma_scaled_delta = (np.max(ap_corr_good) - np.min(ap_corr_good))/psfSigma
898
899 return ap_corr_sigma_scaled_delta
900
901

◆ compute_magnitude_limit()

lsst.pipe.tasks.computeExposureSummaryStats.compute_magnitude_limit ( psfArea,
skyBg,
zeroPoint,
readNoise,
gain,
snr )
Compute the expected point-source magnitude limit at a given
signal-to-noise ratio given the exposure-level metadata. Based on
the signal-to-noise formula provided in SMTN-002 (see LSE-40 for
more details on the calculation).

  SNR = C / sqrt( C/g + (B/g + sigma_inst**2) * neff )

where C is the counts from the source, B is counts from the (sky)
background, sigma_inst is the instrumental (read) noise, neff is
the effective size of the PSF, and g is the gain in e-/ADU.  Note
that input values of ``skyBg``, ``zeroPoint``, and ``readNoise``
should all consistently be in electrons or ADU.

Parameters
----------
psfArea : `float`
    The effective area of the PSF [pix].
skyBg : `float`
    The sky background counts for the exposure [ADU or e-].
zeroPoint : `float`
    The zeropoint (includes exposure time) [ADU or e-].
readNoise : `float`
    The instrumental read noise for the exposure [ADU or e-].
gain : `float`
    The instrumental gain for the exposure [e-/ADU]. The gain should
     be 1.0 if the skyBg, zeroPoint, and readNoise are in e-.
snr : `float`
    Signal-to-noise ratio at which magnitude limit is calculated.

Returns
-------
magnitude_limit : `float`
    The expected magnitude limit at the given signal to noise.

Definition at line 902 of file computeExposureSummaryStats.py.

909):
910 """Compute the expected point-source magnitude limit at a given
911 signal-to-noise ratio given the exposure-level metadata. Based on
912 the signal-to-noise formula provided in SMTN-002 (see LSE-40 for
913 more details on the calculation).
914
915 SNR = C / sqrt( C/g + (B/g + sigma_inst**2) * neff )
916
917 where C is the counts from the source, B is counts from the (sky)
918 background, sigma_inst is the instrumental (read) noise, neff is
919 the effective size of the PSF, and g is the gain in e-/ADU. Note
920 that input values of ``skyBg``, ``zeroPoint``, and ``readNoise``
921 should all consistently be in electrons or ADU.
922
923 Parameters
924 ----------
925 psfArea : `float`
926 The effective area of the PSF [pix].
927 skyBg : `float`
928 The sky background counts for the exposure [ADU or e-].
929 zeroPoint : `float`
930 The zeropoint (includes exposure time) [ADU or e-].
931 readNoise : `float`
932 The instrumental read noise for the exposure [ADU or e-].
933 gain : `float`
934 The instrumental gain for the exposure [e-/ADU]. The gain should
935 be 1.0 if the skyBg, zeroPoint, and readNoise are in e-.
936 snr : `float`
937 Signal-to-noise ratio at which magnitude limit is calculated.
938
939 Returns
940 -------
941 magnitude_limit : `float`
942 The expected magnitude limit at the given signal to noise.
943
944 """
945 # Effective number of pixels within the PSF calculated directly
946 # from the PSF model.
947 neff = psfArea
948
949 # Instrumental noise (read noise only)
950 sigma_inst = readNoise
951
952 # Total background counts derived from Eq. 45 of LSE-40
953 background = (skyBg/gain + sigma_inst**2) * neff
954 # Source counts to achieve the desired SNR (from quadratic formula)
955 source = (snr**2)/(2*gain) + np.sqrt((snr**4)/(4*gain) + snr**2 * background)
956
957 # Convert to a magnitude using the zeropoint.
958 # Note: Zeropoint currently includes exposure time
959 magnitude_limit = -2.5*np.log10(source) + zeroPoint
960
961 return magnitude_limit

◆ compute_psf_image_deltas()

lsst.pipe.tasks.computeExposureSummaryStats.compute_psf_image_deltas ( image_mask,
image_psf,
sampling = 96,
ap_radius_pix = 3.0,
bad_mask_bits = ["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"] )
Compute the delta between the maximum and minimum model PSF trace radius
values evaluated on a grid of points lying in the unmasked region of the
image.

Parameters
----------
image_mask : `lsst.afw.image.Mask`
    The mask plane associated with the exposure.
image_psf : `lsst.afw.detection.Psf`
    The PSF model associated with the exposure.
sampling : `int`, optional
    Sampling rate in each dimension to create the grid of points at which
    to evaluate ``image_psf``s trace radius value. The tradeoff is between
    adequate sampling versus speed.
ap_radius_pix : `float`, optional
    Radius in pixels of the aperture on which to measure the flux of the
    PSF model.
bad_mask_bits : `list` [`str`], optional
    Mask bits required to be absent for a pixel to be considered
    "unmasked".

Returns
-------
psf_trace_radius_delta, psf_ap_flux_delta : `float`
    The delta (in pixels) between the maximum and minimum model PSF trace
    radius values and the PSF aperture fluxes (with aperture radius of
    max(2, 3*psfSigma)) evaluated on the x,y-grid subsampled on the
    unmasked detector pixels by a factor of ``sampling``.  If both the
    model PSF trace radius value and aperture flux value on the grid
    evaluate to NaN, then NaNs are returned immediately.

Definition at line 777 of file computeExposureSummaryStats.py.

783):
784 """Compute the delta between the maximum and minimum model PSF trace radius
785 values evaluated on a grid of points lying in the unmasked region of the
786 image.
787
788 Parameters
789 ----------
790 image_mask : `lsst.afw.image.Mask`
791 The mask plane associated with the exposure.
792 image_psf : `lsst.afw.detection.Psf`
793 The PSF model associated with the exposure.
794 sampling : `int`, optional
795 Sampling rate in each dimension to create the grid of points at which
796 to evaluate ``image_psf``s trace radius value. The tradeoff is between
797 adequate sampling versus speed.
798 ap_radius_pix : `float`, optional
799 Radius in pixels of the aperture on which to measure the flux of the
800 PSF model.
801 bad_mask_bits : `list` [`str`], optional
802 Mask bits required to be absent for a pixel to be considered
803 "unmasked".
804
805 Returns
806 -------
807 psf_trace_radius_delta, psf_ap_flux_delta : `float`
808 The delta (in pixels) between the maximum and minimum model PSF trace
809 radius values and the PSF aperture fluxes (with aperture radius of
810 max(2, 3*psfSigma)) evaluated on the x,y-grid subsampled on the
811 unmasked detector pixels by a factor of ``sampling``. If both the
812 model PSF trace radius value and aperture flux value on the grid
813 evaluate to NaN, then NaNs are returned immediately.
814 """
815 psf_trace_radius_list = []
816 psf_ap_flux_list = []
817 mask_arr = image_mask.array[::sampling, ::sampling]
818 bitmask = image_mask.getPlaneBitMask(bad_mask_bits)
819 good = ((mask_arr & bitmask) == 0)
820
821 x = np.arange(good.shape[1]) * sampling
822 y = np.arange(good.shape[0]) * sampling
823 xx, yy = np.meshgrid(x, y)
824
825 for x_mesh, y_mesh, good_mesh in zip(xx, yy, good):
826 for x_point, y_point, is_good in zip(x_mesh, y_mesh, good_mesh):
827 if is_good:
828 position = geom.Point2D(x_point, y_point)
829 psf_trace_radius = image_psf.computeShape(position).getTraceRadius()
830 psf_ap_flux = image_psf.computeApertureFlux(ap_radius_pix, position)
831 if ~np.isfinite(psf_trace_radius) and ~np.isfinite(psf_ap_flux):
832 return float("nan"), float("nan")
833 psf_trace_radius_list.append(psf_trace_radius)
834 psf_ap_flux_list.append(psf_ap_flux)
835
836 psf_trace_radius_delta = np.max(psf_trace_radius_list) - np.min(psf_trace_radius_list)
837 if np.any(np.asarray(psf_ap_flux_list) < 0.0): # Consider any -ve flux entry as "bad".
838 psf_ap_flux_delta = float("nan")
839 else:
840 psf_ap_flux_delta = np.max(psf_ap_flux_list) - np.min(psf_ap_flux_list)
841
842 return psf_trace_radius_delta, psf_ap_flux_delta
843
844

◆ maximum_nearest_psf_distance()

lsst.pipe.tasks.computeExposureSummaryStats.maximum_nearest_psf_distance ( image_mask,
psf_cat,
sampling = 8,
bad_mask_bits = ["BAD", "CR", "INTRP", "SAT", "SUSPECT", "NO_DATA", "EDGE"] )
Compute the maximum distance of an unmasked pixel to its nearest PSF.

Parameters
----------
image_mask : `lsst.afw.image.Mask`
    The mask plane associated with the exposure.
psf_cat : `lsst.afw.table.SourceCatalog` or `astropy.table.Table`
    Catalog containing only the stars used in the PSF modeling.
sampling : `int`
    Sampling rate in each dimension to create the grid of points on which
    to evaluate the distance to the nearest PSF star. The tradeoff is
    between adequate sampling versus speed.
bad_mask_bits : `list` [`str`]
    Mask bits required to be absent for a pixel to be considered
    "unmasked".

Returns
-------
max_dist_to_nearest_psf : `float`
    The maximum distance (in pixels) of an unmasked pixel to its nearest
    PSF model star.

Definition at line 730 of file computeExposureSummaryStats.py.

735):
736 """Compute the maximum distance of an unmasked pixel to its nearest PSF.
737
738 Parameters
739 ----------
740 image_mask : `lsst.afw.image.Mask`
741 The mask plane associated with the exposure.
742 psf_cat : `lsst.afw.table.SourceCatalog` or `astropy.table.Table`
743 Catalog containing only the stars used in the PSF modeling.
744 sampling : `int`
745 Sampling rate in each dimension to create the grid of points on which
746 to evaluate the distance to the nearest PSF star. The tradeoff is
747 between adequate sampling versus speed.
748 bad_mask_bits : `list` [`str`]
749 Mask bits required to be absent for a pixel to be considered
750 "unmasked".
751
752 Returns
753 -------
754 max_dist_to_nearest_psf : `float`
755 The maximum distance (in pixels) of an unmasked pixel to its nearest
756 PSF model star.
757 """
758 mask_arr = image_mask.array[::sampling, ::sampling]
759 bitmask = image_mask.getPlaneBitMask(bad_mask_bits)
760 good = ((mask_arr & bitmask) == 0)
761
762 x = np.arange(good.shape[1]) * sampling
763 y = np.arange(good.shape[0]) * sampling
764 xx, yy = np.meshgrid(x, y)
765
766 dist_to_nearest_psf = np.full(good.shape, np.inf)
767 for psf in psf_cat:
768 x_psf = psf["slot_Centroid_x"]
769 y_psf = psf["slot_Centroid_y"]
770 dist_to_nearest_psf = np.minimum(dist_to_nearest_psf, np.hypot(xx - x_psf, yy - y_psf))
771 unmasked_dists = dist_to_nearest_psf * good
772 max_dist_to_nearest_psf = np.max(unmasked_dists)
773
774 return max_dist_to_nearest_psf
775
776