LSSTApplications  20.0.0
LSSTDataManagementBasePackage
Functions
lsst.ip.isr.isrFunctions Namespace Reference

Functions

def createPsf (fwhm)
 
def transposeMaskedImage (maskedImage)
 
def interpolateDefectList (maskedImage, defectList, fwhm, fallbackValue=None)
 
def makeThresholdMask (maskedImage, threshold, growFootprints=1, maskName='SAT')
 
def growMasks (mask, radius=0, maskNameList=['BAD'], maskValue="BAD")
 
def interpolateFromMask (maskedImage, fwhm, growSaturatedFootprints=1, maskNameList=['SAT'], fallbackValue=None)
 
def saturationCorrection (maskedImage, saturation, fwhm, growFootprints=1, interpolate=True, maskName='SAT', fallbackValue=None)
 
def trimToMatchCalibBBox (rawMaskedImage, calibMaskedImage)
 
def biasCorrection (maskedImage, biasMaskedImage, trimToFit=False)
 
def darkCorrection (maskedImage, darkMaskedImage, expScale, darkScale, invert=False, trimToFit=False)
 
def updateVariance (maskedImage, gain, readNoise)
 
def flatCorrection (maskedImage, flatMaskedImage, scalingType, userScale=1.0, invert=False, trimToFit=False)
 
def illuminationCorrection (maskedImage, illumMaskedImage, illumScale, trimToFit=True)
 
def overscanCorrection (ampMaskedImage, overscanImage, fitType='MEDIAN', order=1, collapseRej=3.0, statControl=None, overscanIsInt=True)
 
def brighterFatterCorrection (exposure, kernel, maxIter, threshold, applyGain, gains=None)
 
def gainContext (exp, image, apply, gains=None)
 
def attachTransmissionCurve (exposure, opticsTransmission=None, filterTransmission=None, sensorTransmission=None, atmosphereTransmission=None)
 
def applyGains (exposure, normalizeGains=False)
 
def widenSaturationTrails (mask)
 
def setBadRegions (exposure, badStatistic="MEDIAN")
 

Function Documentation

◆ applyGains()

def lsst.ip.isr.isrFunctions.applyGains (   exposure,
  normalizeGains = False 
)
Scale an exposure by the amplifier gains.

Parameters
----------
exposure : `lsst.afw.image.Exposure`
    Exposure to process.  The image is modified.
normalizeGains : `Bool`, optional
    If True, then amplifiers are scaled to force the median of
    each amplifier to equal the median of those medians.

Definition at line 745 of file isrFunctions.py.

745 def applyGains(exposure, normalizeGains=False):
746  """Scale an exposure by the amplifier gains.
747 
748  Parameters
749  ----------
750  exposure : `lsst.afw.image.Exposure`
751  Exposure to process. The image is modified.
752  normalizeGains : `Bool`, optional
753  If True, then amplifiers are scaled to force the median of
754  each amplifier to equal the median of those medians.
755  """
756  ccd = exposure.getDetector()
757  ccdImage = exposure.getMaskedImage()
758 
759  medians = []
760  for amp in ccd:
761  sim = ccdImage.Factory(ccdImage, amp.getBBox())
762  sim *= amp.getGain()
763 
764  if normalizeGains:
765  medians.append(numpy.median(sim.getImage().getArray()))
766 
767  if normalizeGains:
768  median = numpy.median(numpy.array(medians))
769  for index, amp in enumerate(ccd):
770  sim = ccdImage.Factory(ccdImage, amp.getBBox())
771  if medians[index] != 0.0:
772  sim *= median/medians[index]
773 
774 

◆ attachTransmissionCurve()

def lsst.ip.isr.isrFunctions.attachTransmissionCurve (   exposure,
  opticsTransmission = None,
  filterTransmission = None,
  sensorTransmission = None,
  atmosphereTransmission = None 
)
Attach a TransmissionCurve to an Exposure, given separate curves for
different components.

Parameters
----------
exposure : `lsst.afw.image.Exposure`
    Exposure object to modify by attaching the product of all given
    ``TransmissionCurves`` in post-assembly trimmed detector coordinates.
    Must have a valid ``Detector`` attached that matches the detector
    associated with sensorTransmission.
opticsTransmission : `lsst.afw.image.TransmissionCurve`
    A ``TransmissionCurve`` that represents the throughput of the optics,
    to be evaluated in focal-plane coordinates.
filterTransmission : `lsst.afw.image.TransmissionCurve`
    A ``TransmissionCurve`` that represents the throughput of the filter
    itself, to be evaluated in focal-plane coordinates.
sensorTransmission : `lsst.afw.image.TransmissionCurve`
    A ``TransmissionCurve`` that represents the throughput of the sensor
    itself, to be evaluated in post-assembly trimmed detector coordinates.
atmosphereTransmission : `lsst.afw.image.TransmissionCurve`
    A ``TransmissionCurve`` that represents the throughput of the
    atmosphere, assumed to be spatially constant.

Returns
-------
combined : `lsst.afw.image.TransmissionCurve`
    The TransmissionCurve attached to the exposure.

Notes
-----
All ``TransmissionCurve`` arguments are optional; if none are provided, the
attached ``TransmissionCurve`` will have unit transmission everywhere.

Definition at line 693 of file isrFunctions.py.

693 def attachTransmissionCurve(exposure, opticsTransmission=None, filterTransmission=None,
694  sensorTransmission=None, atmosphereTransmission=None):
695  """Attach a TransmissionCurve to an Exposure, given separate curves for
696  different components.
697 
698  Parameters
699  ----------
700  exposure : `lsst.afw.image.Exposure`
701  Exposure object to modify by attaching the product of all given
702  ``TransmissionCurves`` in post-assembly trimmed detector coordinates.
703  Must have a valid ``Detector`` attached that matches the detector
704  associated with sensorTransmission.
705  opticsTransmission : `lsst.afw.image.TransmissionCurve`
706  A ``TransmissionCurve`` that represents the throughput of the optics,
707  to be evaluated in focal-plane coordinates.
708  filterTransmission : `lsst.afw.image.TransmissionCurve`
709  A ``TransmissionCurve`` that represents the throughput of the filter
710  itself, to be evaluated in focal-plane coordinates.
711  sensorTransmission : `lsst.afw.image.TransmissionCurve`
712  A ``TransmissionCurve`` that represents the throughput of the sensor
713  itself, to be evaluated in post-assembly trimmed detector coordinates.
714  atmosphereTransmission : `lsst.afw.image.TransmissionCurve`
715  A ``TransmissionCurve`` that represents the throughput of the
716  atmosphere, assumed to be spatially constant.
717 
718  Returns
719  -------
720  combined : `lsst.afw.image.TransmissionCurve`
721  The TransmissionCurve attached to the exposure.
722 
723  Notes
724  -----
725  All ``TransmissionCurve`` arguments are optional; if none are provided, the
726  attached ``TransmissionCurve`` will have unit transmission everywhere.
727  """
728  combined = afwImage.TransmissionCurve.makeIdentity()
729  if atmosphereTransmission is not None:
730  combined *= atmosphereTransmission
731  if opticsTransmission is not None:
732  combined *= opticsTransmission
733  if filterTransmission is not None:
734  combined *= filterTransmission
735  detector = exposure.getDetector()
736  fpToPix = detector.getTransform(fromSys=camGeom.FOCAL_PLANE,
737  toSys=camGeom.PIXELS)
738  combined = combined.transformedBy(fpToPix)
739  if sensorTransmission is not None:
740  combined *= sensorTransmission
741  exposure.getInfo().setTransmissionCurve(combined)
742  return combined
743 
744 

◆ biasCorrection()

def lsst.ip.isr.isrFunctions.biasCorrection (   maskedImage,
  biasMaskedImage,
  trimToFit = False 
)
Apply bias correction in place.

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
   Image to process.  The image is modified by this method.
biasMaskedImage : `lsst.afw.image.MaskedImage`
    Bias image of the same size as ``maskedImage``
trimToFit : `Bool`, optional
    If True, raw data is symmetrically trimmed to match
    calibration size.

Raises
------
RuntimeError
    Raised if ``maskedImage`` and ``biasMaskedImage`` do not have
    the same size.

Definition at line 269 of file isrFunctions.py.

269 def biasCorrection(maskedImage, biasMaskedImage, trimToFit=False):
270  """Apply bias correction in place.
271 
272  Parameters
273  ----------
274  maskedImage : `lsst.afw.image.MaskedImage`
275  Image to process. The image is modified by this method.
276  biasMaskedImage : `lsst.afw.image.MaskedImage`
277  Bias image of the same size as ``maskedImage``
278  trimToFit : `Bool`, optional
279  If True, raw data is symmetrically trimmed to match
280  calibration size.
281 
282  Raises
283  ------
284  RuntimeError
285  Raised if ``maskedImage`` and ``biasMaskedImage`` do not have
286  the same size.
287 
288  """
289  if trimToFit:
290  maskedImage = trimToMatchCalibBBox(maskedImage, biasMaskedImage)
291 
292  if maskedImage.getBBox(afwImage.LOCAL) != biasMaskedImage.getBBox(afwImage.LOCAL):
293  raise RuntimeError("maskedImage bbox %s != biasMaskedImage bbox %s" %
294  (maskedImage.getBBox(afwImage.LOCAL), biasMaskedImage.getBBox(afwImage.LOCAL)))
295  maskedImage -= biasMaskedImage
296 
297 

◆ brighterFatterCorrection()

def lsst.ip.isr.isrFunctions.brighterFatterCorrection (   exposure,
  kernel,
  maxIter,
  threshold,
  applyGain,
  gains = None 
)
Apply brighter fatter correction in place for the image.

Parameters
----------
exposure : `lsst.afw.image.Exposure`
    Exposure to have brighter-fatter correction applied.  Modified
    by this method.
kernel : `numpy.ndarray`
    Brighter-fatter kernel to apply.
maxIter : scalar
    Number of correction iterations to run.
threshold : scalar
    Convergence threshold in terms of the sum of absolute
    deviations between an iteration and the previous one.
applyGain : `Bool`
    If True, then the exposure values are scaled by the gain prior
    to correction.
gains : `dict` [`str`, `float`]
    A dictionary, keyed by amplifier name, of the gains to use.
    If gains is None, the nominal gains in the amplifier object are used.

Returns
-------
diff : `float`
    Final difference between iterations achieved in correction.
iteration : `int`
    Number of iterations used to calculate correction.

Notes
-----
This correction takes a kernel that has been derived from flat
field images to redistribute the charge.  The gradient of the
kernel is the deflection field due to the accumulated charge.

Given the original image I(x) and the kernel K(x) we can compute
the corrected image Ic(x) using the following equation:

Ic(x) = I(x) + 0.5*d/dx(I(x)*d/dx(int( dy*K(x-y)*I(y))))

To evaluate the derivative term we expand it as follows:

0.5 * ( d/dx(I(x))*d/dx(int(dy*K(x-y)*I(y))) + I(x)*d^2/dx^2(int(dy* K(x-y)*I(y))) )

Because we use the measured counts instead of the incident counts
we apply the correction iteratively to reconstruct the original
counts and the correction.  We stop iterating when the summed
difference between the current corrected image and the one from
the previous iteration is below the threshold.  We do not require
convergence because the number of iterations is too large a
computational cost.  How we define the threshold still needs to be
evaluated, the current default was shown to work reasonably well
on a small set of images.  For more information on the method see
DocuShare Document-19407.

The edges as defined by the kernel are not corrected because they
have spurious values due to the convolution.

Definition at line 518 of file isrFunctions.py.

518 def brighterFatterCorrection(exposure, kernel, maxIter, threshold, applyGain, gains=None):
519  """Apply brighter fatter correction in place for the image.
520 
521  Parameters
522  ----------
523  exposure : `lsst.afw.image.Exposure`
524  Exposure to have brighter-fatter correction applied. Modified
525  by this method.
526  kernel : `numpy.ndarray`
527  Brighter-fatter kernel to apply.
528  maxIter : scalar
529  Number of correction iterations to run.
530  threshold : scalar
531  Convergence threshold in terms of the sum of absolute
532  deviations between an iteration and the previous one.
533  applyGain : `Bool`
534  If True, then the exposure values are scaled by the gain prior
535  to correction.
536  gains : `dict` [`str`, `float`]
537  A dictionary, keyed by amplifier name, of the gains to use.
538  If gains is None, the nominal gains in the amplifier object are used.
539 
540  Returns
541  -------
542  diff : `float`
543  Final difference between iterations achieved in correction.
544  iteration : `int`
545  Number of iterations used to calculate correction.
546 
547  Notes
548  -----
549  This correction takes a kernel that has been derived from flat
550  field images to redistribute the charge. The gradient of the
551  kernel is the deflection field due to the accumulated charge.
552 
553  Given the original image I(x) and the kernel K(x) we can compute
554  the corrected image Ic(x) using the following equation:
555 
556  Ic(x) = I(x) + 0.5*d/dx(I(x)*d/dx(int( dy*K(x-y)*I(y))))
557 
558  To evaluate the derivative term we expand it as follows:
559 
560  0.5 * ( d/dx(I(x))*d/dx(int(dy*K(x-y)*I(y))) + I(x)*d^2/dx^2(int(dy* K(x-y)*I(y))) )
561 
562  Because we use the measured counts instead of the incident counts
563  we apply the correction iteratively to reconstruct the original
564  counts and the correction. We stop iterating when the summed
565  difference between the current corrected image and the one from
566  the previous iteration is below the threshold. We do not require
567  convergence because the number of iterations is too large a
568  computational cost. How we define the threshold still needs to be
569  evaluated, the current default was shown to work reasonably well
570  on a small set of images. For more information on the method see
571  DocuShare Document-19407.
572 
573  The edges as defined by the kernel are not corrected because they
574  have spurious values due to the convolution.
575  """
576  image = exposure.getMaskedImage().getImage()
577 
578  # The image needs to be units of electrons/holes
579  with gainContext(exposure, image, applyGain, gains):
580 
581  kLx = numpy.shape(kernel)[0]
582  kLy = numpy.shape(kernel)[1]
583  kernelImage = afwImage.ImageD(kLx, kLy)
584  kernelImage.getArray()[:, :] = kernel
585  tempImage = image.clone()
586 
587  nanIndex = numpy.isnan(tempImage.getArray())
588  tempImage.getArray()[nanIndex] = 0.
589 
590  outImage = afwImage.ImageF(image.getDimensions())
591  corr = numpy.zeros_like(image.getArray())
592  prev_image = numpy.zeros_like(image.getArray())
593  convCntrl = afwMath.ConvolutionControl(False, True, 1)
594  fixedKernel = afwMath.FixedKernel(kernelImage)
595 
596  # Define boundary by convolution region. The region that the correction will be
597  # calculated for is one fewer in each dimension because of the second derivative terms.
598  # NOTE: these need to use integer math, as we're using start:end as numpy index ranges.
599  startX = kLx//2
600  endX = -kLx//2
601  startY = kLy//2
602  endY = -kLy//2
603 
604  for iteration in range(maxIter):
605 
606  afwMath.convolve(outImage, tempImage, fixedKernel, convCntrl)
607  tmpArray = tempImage.getArray()
608  outArray = outImage.getArray()
609 
610  with numpy.errstate(invalid="ignore", over="ignore"):
611  # First derivative term
612  gradTmp = numpy.gradient(tmpArray[startY:endY, startX:endX])
613  gradOut = numpy.gradient(outArray[startY:endY, startX:endX])
614  first = (gradTmp[0]*gradOut[0] + gradTmp[1]*gradOut[1])[1:-1, 1:-1]
615 
616  # Second derivative term
617  diffOut20 = numpy.diff(outArray, 2, 0)[startY:endY, startX + 1:endX - 1]
618  diffOut21 = numpy.diff(outArray, 2, 1)[startY + 1:endY - 1, startX:endX]
619  second = tmpArray[startY + 1:endY - 1, startX + 1:endX - 1]*(diffOut20 + diffOut21)
620 
621  corr[startY + 1:endY - 1, startX + 1:endX - 1] = 0.5*(first + second)
622 
623  tmpArray[:, :] = image.getArray()[:, :]
624  tmpArray[nanIndex] = 0.
625  tmpArray[startY:endY, startX:endX] += corr[startY:endY, startX:endX]
626 
627  if iteration > 0:
628  diff = numpy.sum(numpy.abs(prev_image - tmpArray))
629 
630  if diff < threshold:
631  break
632  prev_image[:, :] = tmpArray[:, :]
633 
634  image.getArray()[startY + 1:endY - 1, startX + 1:endX - 1] += \
635  corr[startY + 1:endY - 1, startX + 1:endX - 1]
636 
637  return diff, iteration
638 
639 
640 @contextmanager

◆ createPsf()

def lsst.ip.isr.isrFunctions.createPsf (   fwhm)
Make a double Gaussian PSF.

Parameters
----------
fwhm : scalar
    FWHM of double Gaussian smoothing kernel.

Returns
-------
psf : `lsst.meas.algorithms.DoubleGaussianPsf`
    The created smoothing kernel.

Definition at line 39 of file isrFunctions.py.

39 def createPsf(fwhm):
40  """Make a double Gaussian PSF.
41 
42  Parameters
43  ----------
44  fwhm : scalar
45  FWHM of double Gaussian smoothing kernel.
46 
47  Returns
48  -------
49  psf : `lsst.meas.algorithms.DoubleGaussianPsf`
50  The created smoothing kernel.
51  """
52  ksize = 4*int(fwhm) + 1
53  return measAlg.DoubleGaussianPsf(ksize, ksize, fwhm/(2*math.sqrt(2*math.log(2))))
54 
55 

◆ darkCorrection()

def lsst.ip.isr.isrFunctions.darkCorrection (   maskedImage,
  darkMaskedImage,
  expScale,
  darkScale,
  invert = False,
  trimToFit = False 
)
Apply dark correction in place.

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
   Image to process.  The image is modified by this method.
darkMaskedImage : `lsst.afw.image.MaskedImage`
    Dark image of the same size as ``maskedImage``.
expScale : scalar
    Dark exposure time for ``maskedImage``.
darkScale : scalar
    Dark exposure time for ``darkMaskedImage``.
invert : `Bool`, optional
    If True, re-add the dark to an already corrected image.
trimToFit : `Bool`, optional
    If True, raw data is symmetrically trimmed to match
    calibration size.

Raises
------
RuntimeError
    Raised if ``maskedImage`` and ``darkMaskedImage`` do not have
    the same size.

Notes
-----
The dark correction is applied by calculating:
    maskedImage -= dark * expScaling / darkScaling

Definition at line 298 of file isrFunctions.py.

298 def darkCorrection(maskedImage, darkMaskedImage, expScale, darkScale, invert=False, trimToFit=False):
299  """Apply dark correction in place.
300 
301  Parameters
302  ----------
303  maskedImage : `lsst.afw.image.MaskedImage`
304  Image to process. The image is modified by this method.
305  darkMaskedImage : `lsst.afw.image.MaskedImage`
306  Dark image of the same size as ``maskedImage``.
307  expScale : scalar
308  Dark exposure time for ``maskedImage``.
309  darkScale : scalar
310  Dark exposure time for ``darkMaskedImage``.
311  invert : `Bool`, optional
312  If True, re-add the dark to an already corrected image.
313  trimToFit : `Bool`, optional
314  If True, raw data is symmetrically trimmed to match
315  calibration size.
316 
317  Raises
318  ------
319  RuntimeError
320  Raised if ``maskedImage`` and ``darkMaskedImage`` do not have
321  the same size.
322 
323  Notes
324  -----
325  The dark correction is applied by calculating:
326  maskedImage -= dark * expScaling / darkScaling
327  """
328  if trimToFit:
329  maskedImage = trimToMatchCalibBBox(maskedImage, darkMaskedImage)
330 
331  if maskedImage.getBBox(afwImage.LOCAL) != darkMaskedImage.getBBox(afwImage.LOCAL):
332  raise RuntimeError("maskedImage bbox %s != darkMaskedImage bbox %s" %
333  (maskedImage.getBBox(afwImage.LOCAL), darkMaskedImage.getBBox(afwImage.LOCAL)))
334 
335  scale = expScale / darkScale
336  if not invert:
337  maskedImage.scaledMinus(scale, darkMaskedImage)
338  else:
339  maskedImage.scaledPlus(scale, darkMaskedImage)
340 
341 

◆ flatCorrection()

def lsst.ip.isr.isrFunctions.flatCorrection (   maskedImage,
  flatMaskedImage,
  scalingType,
  userScale = 1.0,
  invert = False,
  trimToFit = False 
)
Apply flat correction in place.

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
    Image to process.  The image is modified.
flatMaskedImage : `lsst.afw.image.MaskedImage`
    Flat image of the same size as ``maskedImage``
scalingType : str
    Flat scale computation method.  Allowed values are 'MEAN',
    'MEDIAN', or 'USER'.
userScale : scalar, optional
    Scale to use if ``scalingType``='USER'.
invert : `Bool`, optional
    If True, unflatten an already flattened image.
trimToFit : `Bool`, optional
    If True, raw data is symmetrically trimmed to match
    calibration size.

Raises
------
RuntimeError
    Raised if ``maskedImage`` and ``flatMaskedImage`` do not have
    the same size or if ``scalingType`` is not an allowed value.

Definition at line 360 of file isrFunctions.py.

360 def flatCorrection(maskedImage, flatMaskedImage, scalingType, userScale=1.0, invert=False, trimToFit=False):
361  """Apply flat correction in place.
362 
363  Parameters
364  ----------
365  maskedImage : `lsst.afw.image.MaskedImage`
366  Image to process. The image is modified.
367  flatMaskedImage : `lsst.afw.image.MaskedImage`
368  Flat image of the same size as ``maskedImage``
369  scalingType : str
370  Flat scale computation method. Allowed values are 'MEAN',
371  'MEDIAN', or 'USER'.
372  userScale : scalar, optional
373  Scale to use if ``scalingType``='USER'.
374  invert : `Bool`, optional
375  If True, unflatten an already flattened image.
376  trimToFit : `Bool`, optional
377  If True, raw data is symmetrically trimmed to match
378  calibration size.
379 
380  Raises
381  ------
382  RuntimeError
383  Raised if ``maskedImage`` and ``flatMaskedImage`` do not have
384  the same size or if ``scalingType`` is not an allowed value.
385  """
386  if trimToFit:
387  maskedImage = trimToMatchCalibBBox(maskedImage, flatMaskedImage)
388 
389  if maskedImage.getBBox(afwImage.LOCAL) != flatMaskedImage.getBBox(afwImage.LOCAL):
390  raise RuntimeError("maskedImage bbox %s != flatMaskedImage bbox %s" %
391  (maskedImage.getBBox(afwImage.LOCAL), flatMaskedImage.getBBox(afwImage.LOCAL)))
392 
393  # Figure out scale from the data
394  # Ideally the flats are normalized by the calibration product pipeline, but this allows some flexibility
395  # in the case that the flat is created by some other mechanism.
396  if scalingType in ('MEAN', 'MEDIAN'):
397  scalingType = afwMath.stringToStatisticsProperty(scalingType)
398  flatScale = afwMath.makeStatistics(flatMaskedImage.image, scalingType).getValue()
399  elif scalingType == 'USER':
400  flatScale = userScale
401  else:
402  raise RuntimeError('%s : %s not implemented' % ("flatCorrection", scalingType))
403 
404  if not invert:
405  maskedImage.scaledDivides(1.0/flatScale, flatMaskedImage)
406  else:
407  maskedImage.scaledMultiplies(1.0/flatScale, flatMaskedImage)
408 
409 

◆ gainContext()

def lsst.ip.isr.isrFunctions.gainContext (   exp,
  image,
  apply,
  gains = None 
)
Context manager that applies and removes gain.

Parameters
----------
exp : `lsst.afw.image.Exposure`
    Exposure to apply/remove gain.
image : `lsst.afw.image.Image`
    Image to apply/remove gain.
apply : `Bool`
    If True, apply and remove the amplifier gain.
gains : `dict` [`str`, `float`]
    A dictionary, keyed by amplifier name, of the gains to use.
    If gains is None, the nominal gains in the amplifier object are used.

Yields
------
exp : `lsst.afw.image.Exposure`
    Exposure with the gain applied.

Definition at line 641 of file isrFunctions.py.

641 def gainContext(exp, image, apply, gains=None):
642  """Context manager that applies and removes gain.
643 
644  Parameters
645  ----------
646  exp : `lsst.afw.image.Exposure`
647  Exposure to apply/remove gain.
648  image : `lsst.afw.image.Image`
649  Image to apply/remove gain.
650  apply : `Bool`
651  If True, apply and remove the amplifier gain.
652  gains : `dict` [`str`, `float`]
653  A dictionary, keyed by amplifier name, of the gains to use.
654  If gains is None, the nominal gains in the amplifier object are used.
655 
656  Yields
657  ------
658  exp : `lsst.afw.image.Exposure`
659  Exposure with the gain applied.
660  """
661  # check we have all of them if provided because mixing and matching would
662  # be a real mess
663  if gains and apply is True:
664  ampNames = [amp.getName() for amp in exp.getDetector()]
665  for ampName in ampNames:
666  if ampName not in gains.keys():
667  raise RuntimeError(f"Gains provided to gain context, but no entry found for amp {ampName}")
668 
669  if apply:
670  ccd = exp.getDetector()
671  for amp in ccd:
672  sim = image.Factory(image, amp.getBBox())
673  if gains:
674  gain = gains[amp.getName()]
675  else:
676  gain = amp.getGain()
677  sim *= gain
678 
679  try:
680  yield exp
681  finally:
682  if apply:
683  ccd = exp.getDetector()
684  for amp in ccd:
685  sim = image.Factory(image, amp.getBBox())
686  if gains:
687  gain = gains[amp.getName()]
688  else:
689  gain = amp.getGain()
690  sim /= gain
691 
692 

◆ growMasks()

def lsst.ip.isr.isrFunctions.growMasks (   mask,
  radius = 0,
  maskNameList = ['BAD'],
  maskValue = "BAD" 
)
Grow a mask by an amount and add to the requested plane.

Parameters
----------
mask : `lsst.afw.image.Mask`
    Mask image to process.
radius : scalar
    Amount to grow the mask.
maskNameList : `str` or `list` [`str`]
    Mask names that should be grown.
maskValue : `str`
    Mask plane to assign the newly masked pixels to.

Definition at line 135 of file isrFunctions.py.

135 def growMasks(mask, radius=0, maskNameList=['BAD'], maskValue="BAD"):
136  """Grow a mask by an amount and add to the requested plane.
137 
138  Parameters
139  ----------
140  mask : `lsst.afw.image.Mask`
141  Mask image to process.
142  radius : scalar
143  Amount to grow the mask.
144  maskNameList : `str` or `list` [`str`]
145  Mask names that should be grown.
146  maskValue : `str`
147  Mask plane to assign the newly masked pixels to.
148  """
149  if radius > 0:
150  thresh = afwDetection.Threshold(mask.getPlaneBitMask(maskNameList), afwDetection.Threshold.BITMASK)
151  fpSet = afwDetection.FootprintSet(mask, thresh)
152  fpSet = afwDetection.FootprintSet(fpSet, rGrow=radius, isotropic=False)
153  fpSet.setMask(mask, maskValue)
154 
155 

◆ illuminationCorrection()

def lsst.ip.isr.isrFunctions.illuminationCorrection (   maskedImage,
  illumMaskedImage,
  illumScale,
  trimToFit = True 
)
Apply illumination correction in place.

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
    Image to process.  The image is modified.
illumMaskedImage : `lsst.afw.image.MaskedImage`
    Illumination correction image of the same size as ``maskedImage``.
illumScale : scalar
    Scale factor for the illumination correction.
trimToFit : `Bool`, optional
    If True, raw data is symmetrically trimmed to match
    calibration size.

Raises
------
RuntimeError
    Raised if ``maskedImage`` and ``illumMaskedImage`` do not have
    the same size.

Definition at line 410 of file isrFunctions.py.

410 def illuminationCorrection(maskedImage, illumMaskedImage, illumScale, trimToFit=True):
411  """Apply illumination correction in place.
412 
413  Parameters
414  ----------
415  maskedImage : `lsst.afw.image.MaskedImage`
416  Image to process. The image is modified.
417  illumMaskedImage : `lsst.afw.image.MaskedImage`
418  Illumination correction image of the same size as ``maskedImage``.
419  illumScale : scalar
420  Scale factor for the illumination correction.
421  trimToFit : `Bool`, optional
422  If True, raw data is symmetrically trimmed to match
423  calibration size.
424 
425  Raises
426  ------
427  RuntimeError
428  Raised if ``maskedImage`` and ``illumMaskedImage`` do not have
429  the same size.
430  """
431  if trimToFit:
432  maskedImage = trimToMatchCalibBBox(maskedImage, illumMaskedImage)
433 
434  if maskedImage.getBBox(afwImage.LOCAL) != illumMaskedImage.getBBox(afwImage.LOCAL):
435  raise RuntimeError("maskedImage bbox %s != illumMaskedImage bbox %s" %
436  (maskedImage.getBBox(afwImage.LOCAL), illumMaskedImage.getBBox(afwImage.LOCAL)))
437 
438  maskedImage.scaledDivides(1.0/illumScale, illumMaskedImage)
439 
440 

◆ interpolateDefectList()

def lsst.ip.isr.isrFunctions.interpolateDefectList (   maskedImage,
  defectList,
  fwhm,
  fallbackValue = None 
)
Interpolate over defects specified in a defect list.

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
    Image to process.
defectList : `lsst.meas.algorithms.Defects`
    List of defects to interpolate over.
fwhm : scalar
    FWHM of double Gaussian smoothing kernel.
fallbackValue : scalar, optional
    Fallback value if an interpolated value cannot be determined.
    If None, then the clipped mean of the image is used.

Definition at line 76 of file isrFunctions.py.

76 def interpolateDefectList(maskedImage, defectList, fwhm, fallbackValue=None):
77  """Interpolate over defects specified in a defect list.
78 
79  Parameters
80  ----------
81  maskedImage : `lsst.afw.image.MaskedImage`
82  Image to process.
83  defectList : `lsst.meas.algorithms.Defects`
84  List of defects to interpolate over.
85  fwhm : scalar
86  FWHM of double Gaussian smoothing kernel.
87  fallbackValue : scalar, optional
88  Fallback value if an interpolated value cannot be determined.
89  If None, then the clipped mean of the image is used.
90  """
91  psf = createPsf(fwhm)
92  if fallbackValue is None:
93  fallbackValue = afwMath.makeStatistics(maskedImage.getImage(), afwMath.MEANCLIP).getValue()
94  if 'INTRP' not in maskedImage.getMask().getMaskPlaneDict():
95  maskedImage.getMask().addMaskPlane('INTRP')
96  measAlg.interpolateOverDefects(maskedImage, psf, defectList, fallbackValue, True)
97  return maskedImage
98 
99 

◆ interpolateFromMask()

def lsst.ip.isr.isrFunctions.interpolateFromMask (   maskedImage,
  fwhm,
  growSaturatedFootprints = 1,
  maskNameList = ['SAT'],
  fallbackValue = None 
)
Interpolate over defects identified by a particular set of mask planes.

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
    Image to process.
fwhm : scalar
    FWHM of double Gaussian smoothing kernel.
growSaturatedFootprints : scalar, optional
    Number of pixels to grow footprints for saturated pixels.
maskNameList : `List` of `str`, optional
    Mask plane name.
fallbackValue : scalar, optional
    Value of last resort for interpolation.

Definition at line 156 of file isrFunctions.py.

156 def interpolateFromMask(maskedImage, fwhm, growSaturatedFootprints=1,
157  maskNameList=['SAT'], fallbackValue=None):
158  """Interpolate over defects identified by a particular set of mask planes.
159 
160  Parameters
161  ----------
162  maskedImage : `lsst.afw.image.MaskedImage`
163  Image to process.
164  fwhm : scalar
165  FWHM of double Gaussian smoothing kernel.
166  growSaturatedFootprints : scalar, optional
167  Number of pixels to grow footprints for saturated pixels.
168  maskNameList : `List` of `str`, optional
169  Mask plane name.
170  fallbackValue : scalar, optional
171  Value of last resort for interpolation.
172  """
173  mask = maskedImage.getMask()
174 
175  if growSaturatedFootprints > 0 and "SAT" in maskNameList:
176  # If we are interpolating over an area larger than the original masked region, we need
177  # to expand the original mask bit to the full area to explain why we interpolated there.
178  growMasks(mask, radius=growSaturatedFootprints, maskNameList=['SAT'], maskValue="SAT")
179 
180  thresh = afwDetection.Threshold(mask.getPlaneBitMask(maskNameList), afwDetection.Threshold.BITMASK)
181  fpSet = afwDetection.FootprintSet(mask, thresh)
182  defectList = measAlg.Defects.fromFootprintList(fpSet.getFootprints())
183 
184  interpolateDefectList(maskedImage, defectList, fwhm, fallbackValue=fallbackValue)
185 
186  return maskedImage
187 
188 

◆ makeThresholdMask()

def lsst.ip.isr.isrFunctions.makeThresholdMask (   maskedImage,
  threshold,
  growFootprints = 1,
  maskName = 'SAT' 
)
Mask pixels based on threshold detection.

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
    Image to process.  Only the mask plane is updated.
threshold : scalar
    Detection threshold.
growFootprints : scalar, optional
    Number of pixels to grow footprints of detected regions.
maskName : str, optional
    Mask plane name, or list of names to convert

Returns
-------
defectList : `lsst.meas.algorithms.Defects`
    Defect list constructed from pixels above the threshold.

Definition at line 100 of file isrFunctions.py.

100 def makeThresholdMask(maskedImage, threshold, growFootprints=1, maskName='SAT'):
101  """Mask pixels based on threshold detection.
102 
103  Parameters
104  ----------
105  maskedImage : `lsst.afw.image.MaskedImage`
106  Image to process. Only the mask plane is updated.
107  threshold : scalar
108  Detection threshold.
109  growFootprints : scalar, optional
110  Number of pixels to grow footprints of detected regions.
111  maskName : str, optional
112  Mask plane name, or list of names to convert
113 
114  Returns
115  -------
116  defectList : `lsst.meas.algorithms.Defects`
117  Defect list constructed from pixels above the threshold.
118  """
119  # find saturated regions
120  thresh = afwDetection.Threshold(threshold)
121  fs = afwDetection.FootprintSet(maskedImage, thresh)
122 
123  if growFootprints > 0:
124  fs = afwDetection.FootprintSet(fs, rGrow=growFootprints, isotropic=False)
125  fpList = fs.getFootprints()
126 
127  # set mask
128  mask = maskedImage.getMask()
129  bitmask = mask.getPlaneBitMask(maskName)
130  afwDetection.setMaskFromFootprintList(mask, fpList, bitmask)
131 
132  return measAlg.Defects.fromFootprintList(fpList)
133 
134 

◆ overscanCorrection()

def lsst.ip.isr.isrFunctions.overscanCorrection (   ampMaskedImage,
  overscanImage,
  fitType = 'MEDIAN',
  order = 1,
  collapseRej = 3.0,
  statControl = None,
  overscanIsInt = True 
)
Apply overscan correction in place.

Parameters
----------
ampMaskedImage : `lsst.afw.image.MaskedImage`
    Image of amplifier to correct; modified.
overscanImage : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
    Image of overscan; modified.
fitType : `str`
    Type of fit for overscan correction. May be one of:

    - ``MEAN``: use mean of overscan.
    - ``MEANCLIP``: use clipped mean of overscan.
    - ``MEDIAN``: use median of overscan.
    - ``MEDIAN_PER_ROW``: use median per row of overscan.
    - ``POLY``: fit with ordinary polynomial.
    - ``CHEB``: fit with Chebyshev polynomial.
    - ``LEG``: fit with Legendre polynomial.
    - ``NATURAL_SPLINE``: fit with natural spline.
    - ``CUBIC_SPLINE``: fit with cubic spline.
    - ``AKIMA_SPLINE``: fit with Akima spline.

order : `int`
    Polynomial order or number of spline knots; ignored unless
    ``fitType`` indicates a polynomial or spline.
statControl : `lsst.afw.math.StatisticsControl`
    Statistics control object.  In particular, we pay attention to numSigmaClip
overscanIsInt : `bool`
    Treat the overscan region as consisting of integers, even if it's been
    converted to float.  E.g. handle ties properly.

Returns
-------
result : `lsst.pipe.base.Struct`
    Result struct with components:

    - ``imageFit``: Value(s) removed from image (scalar or
        `lsst.afw.image.Image`)
    - ``overscanFit``: Value(s) removed from overscan (scalar or
        `lsst.afw.image.Image`)
    - ``overscanImage``: Overscan corrected overscan region
        (`lsst.afw.image.Image`)
Raises
------
pexExcept.Exception
    Raised if ``fitType`` is not an allowed value.

Notes
-----
The ``ampMaskedImage`` and ``overscanImage`` are modified, with the fit
subtracted. Note that the ``overscanImage`` should not be a subimage of
the ``ampMaskedImage``, to avoid being subtracted twice.

Debug plots are available for the SPLINE fitTypes by setting the
`debug.display` for `name` == "lsst.ip.isr.isrFunctions".  These
plots show the scatter plot of the overscan data (collapsed along
the perpendicular dimension) as a function of position on the CCD
(normalized between +/-1).

Definition at line 441 of file isrFunctions.py.

441 def overscanCorrection(ampMaskedImage, overscanImage, fitType='MEDIAN', order=1, collapseRej=3.0,
442  statControl=None, overscanIsInt=True):
443  """Apply overscan correction in place.
444 
445  Parameters
446  ----------
447  ampMaskedImage : `lsst.afw.image.MaskedImage`
448  Image of amplifier to correct; modified.
449  overscanImage : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage`
450  Image of overscan; modified.
451  fitType : `str`
452  Type of fit for overscan correction. May be one of:
453 
454  - ``MEAN``: use mean of overscan.
455  - ``MEANCLIP``: use clipped mean of overscan.
456  - ``MEDIAN``: use median of overscan.
457  - ``MEDIAN_PER_ROW``: use median per row of overscan.
458  - ``POLY``: fit with ordinary polynomial.
459  - ``CHEB``: fit with Chebyshev polynomial.
460  - ``LEG``: fit with Legendre polynomial.
461  - ``NATURAL_SPLINE``: fit with natural spline.
462  - ``CUBIC_SPLINE``: fit with cubic spline.
463  - ``AKIMA_SPLINE``: fit with Akima spline.
464 
465  order : `int`
466  Polynomial order or number of spline knots; ignored unless
467  ``fitType`` indicates a polynomial or spline.
468  statControl : `lsst.afw.math.StatisticsControl`
469  Statistics control object. In particular, we pay attention to numSigmaClip
470  overscanIsInt : `bool`
471  Treat the overscan region as consisting of integers, even if it's been
472  converted to float. E.g. handle ties properly.
473 
474  Returns
475  -------
476  result : `lsst.pipe.base.Struct`
477  Result struct with components:
478 
479  - ``imageFit``: Value(s) removed from image (scalar or
480  `lsst.afw.image.Image`)
481  - ``overscanFit``: Value(s) removed from overscan (scalar or
482  `lsst.afw.image.Image`)
483  - ``overscanImage``: Overscan corrected overscan region
484  (`lsst.afw.image.Image`)
485  Raises
486  ------
487  pexExcept.Exception
488  Raised if ``fitType`` is not an allowed value.
489 
490  Notes
491  -----
492  The ``ampMaskedImage`` and ``overscanImage`` are modified, with the fit
493  subtracted. Note that the ``overscanImage`` should not be a subimage of
494  the ``ampMaskedImage``, to avoid being subtracted twice.
495 
496  Debug plots are available for the SPLINE fitTypes by setting the
497  `debug.display` for `name` == "lsst.ip.isr.isrFunctions". These
498  plots show the scatter plot of the overscan data (collapsed along
499  the perpendicular dimension) as a function of position on the CCD
500  (normalized between +/-1).
501  """
502  ampImage = ampMaskedImage.getImage()
503 
504  config = OverscanCorrectionTaskConfig()
505  if fitType:
506  config.fitType = fitType
507  if order:
508  config.order = order
509  if collapseRej:
510  config.numSigmaClip = collapseRej
511  if overscanIsInt:
512  config.overscanIsInt = True
513 
514  overscanTask = OverscanCorrectionTask(config=config)
515  return overscanTask.run(ampImage, overscanImage)
516 
517 

◆ saturationCorrection()

def lsst.ip.isr.isrFunctions.saturationCorrection (   maskedImage,
  saturation,
  fwhm,
  growFootprints = 1,
  interpolate = True,
  maskName = 'SAT',
  fallbackValue = None 
)
Mark saturated pixels and optionally interpolate over them

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
    Image to process.
saturation  : scalar
    Saturation level used as the detection threshold.
fwhm : scalar
    FWHM of double Gaussian smoothing kernel.
growFootprints : scalar, optional
    Number of pixels to grow footprints of detected regions.
interpolate : Bool, optional
    If True, saturated pixels are interpolated over.
maskName : str, optional
    Mask plane name.
fallbackValue : scalar, optional
    Value of last resort for interpolation.

Definition at line 189 of file isrFunctions.py.

189 def saturationCorrection(maskedImage, saturation, fwhm, growFootprints=1, interpolate=True, maskName='SAT',
190  fallbackValue=None):
191  """Mark saturated pixels and optionally interpolate over them
192 
193  Parameters
194  ----------
195  maskedImage : `lsst.afw.image.MaskedImage`
196  Image to process.
197  saturation : scalar
198  Saturation level used as the detection threshold.
199  fwhm : scalar
200  FWHM of double Gaussian smoothing kernel.
201  growFootprints : scalar, optional
202  Number of pixels to grow footprints of detected regions.
203  interpolate : Bool, optional
204  If True, saturated pixels are interpolated over.
205  maskName : str, optional
206  Mask plane name.
207  fallbackValue : scalar, optional
208  Value of last resort for interpolation.
209  """
210  defectList = makeThresholdMask(
211  maskedImage=maskedImage,
212  threshold=saturation,
213  growFootprints=growFootprints,
214  maskName=maskName,
215  )
216  if interpolate:
217  interpolateDefectList(maskedImage, defectList, fwhm, fallbackValue=fallbackValue)
218 
219  return maskedImage
220 
221 

◆ setBadRegions()

def lsst.ip.isr.isrFunctions.setBadRegions (   exposure,
  badStatistic = "MEDIAN" 
)
Set all BAD areas of the chip to the average of the rest of the exposure

Parameters
----------
exposure : `lsst.afw.image.Exposure`
    Exposure to mask.  The exposure mask is modified.
badStatistic : `str`, optional
    Statistic to use to generate the replacement value from the
    image data.  Allowed values are 'MEDIAN' or 'MEANCLIP'.

Returns
-------
badPixelCount : scalar
    Number of bad pixels masked.
badPixelValue : scalar
    Value substituted for bad pixels.

Raises
------
RuntimeError
    Raised if `badStatistic` is not an allowed value.

Definition at line 822 of file isrFunctions.py.

822 def setBadRegions(exposure, badStatistic="MEDIAN"):
823  """Set all BAD areas of the chip to the average of the rest of the exposure
824 
825  Parameters
826  ----------
827  exposure : `lsst.afw.image.Exposure`
828  Exposure to mask. The exposure mask is modified.
829  badStatistic : `str`, optional
830  Statistic to use to generate the replacement value from the
831  image data. Allowed values are 'MEDIAN' or 'MEANCLIP'.
832 
833  Returns
834  -------
835  badPixelCount : scalar
836  Number of bad pixels masked.
837  badPixelValue : scalar
838  Value substituted for bad pixels.
839 
840  Raises
841  ------
842  RuntimeError
843  Raised if `badStatistic` is not an allowed value.
844  """
845  if badStatistic == "MEDIAN":
846  statistic = afwMath.MEDIAN
847  elif badStatistic == "MEANCLIP":
848  statistic = afwMath.MEANCLIP
849  else:
850  raise RuntimeError("Impossible method %s of bad region correction" % badStatistic)
851 
852  mi = exposure.getMaskedImage()
853  mask = mi.getMask()
854  BAD = mask.getPlaneBitMask("BAD")
855  INTRP = mask.getPlaneBitMask("INTRP")
856 
857  sctrl = afwMath.StatisticsControl()
858  sctrl.setAndMask(BAD)
859  value = afwMath.makeStatistics(mi, statistic, sctrl).getValue()
860 
861  maskArray = mask.getArray()
862  imageArray = mi.getImage().getArray()
863  badPixels = numpy.logical_and((maskArray & BAD) > 0, (maskArray & INTRP) == 0)
864  imageArray[:] = numpy.where(badPixels, value, imageArray)
865 
866  return badPixels.sum(), value

◆ transposeMaskedImage()

def lsst.ip.isr.isrFunctions.transposeMaskedImage (   maskedImage)
Make a transposed copy of a masked image.

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
    Image to process.

Returns
-------
transposed : `lsst.afw.image.MaskedImage`
    The transposed copy of the input image.

Definition at line 56 of file isrFunctions.py.

56 def transposeMaskedImage(maskedImage):
57  """Make a transposed copy of a masked image.
58 
59  Parameters
60  ----------
61  maskedImage : `lsst.afw.image.MaskedImage`
62  Image to process.
63 
64  Returns
65  -------
66  transposed : `lsst.afw.image.MaskedImage`
67  The transposed copy of the input image.
68  """
69  transposed = maskedImage.Factory(lsst.geom.Extent2I(maskedImage.getHeight(), maskedImage.getWidth()))
70  transposed.getImage().getArray()[:] = maskedImage.getImage().getArray().T
71  transposed.getMask().getArray()[:] = maskedImage.getMask().getArray().T
72  transposed.getVariance().getArray()[:] = maskedImage.getVariance().getArray().T
73  return transposed
74 
75 

◆ trimToMatchCalibBBox()

def lsst.ip.isr.isrFunctions.trimToMatchCalibBBox (   rawMaskedImage,
  calibMaskedImage 
)
Compute number of edge trim pixels to match the calibration data.

Use the dimension difference between the raw exposure and the
calibration exposure to compute the edge trim pixels.  This trim
is applied symmetrically, with the same number of pixels masked on
each side.

Parameters
----------
rawMaskedImage : `lsst.afw.image.MaskedImage`
    Image to trim.
calibMaskedImage : `lsst.afw.image.MaskedImage`
    Calibration image to draw new bounding box from.

Returns
-------
replacementMaskedImage : `lsst.afw.image.MaskedImage`
    ``rawMaskedImage`` trimmed to the appropriate size
Raises
------
RuntimeError
   Rasied if ``rawMaskedImage`` cannot be symmetrically trimmed to
   match ``calibMaskedImage``.

Definition at line 222 of file isrFunctions.py.

222 def trimToMatchCalibBBox(rawMaskedImage, calibMaskedImage):
223  """Compute number of edge trim pixels to match the calibration data.
224 
225  Use the dimension difference between the raw exposure and the
226  calibration exposure to compute the edge trim pixels. This trim
227  is applied symmetrically, with the same number of pixels masked on
228  each side.
229 
230  Parameters
231  ----------
232  rawMaskedImage : `lsst.afw.image.MaskedImage`
233  Image to trim.
234  calibMaskedImage : `lsst.afw.image.MaskedImage`
235  Calibration image to draw new bounding box from.
236 
237  Returns
238  -------
239  replacementMaskedImage : `lsst.afw.image.MaskedImage`
240  ``rawMaskedImage`` trimmed to the appropriate size
241  Raises
242  ------
243  RuntimeError
244  Rasied if ``rawMaskedImage`` cannot be symmetrically trimmed to
245  match ``calibMaskedImage``.
246  """
247  nx, ny = rawMaskedImage.getBBox().getDimensions() - calibMaskedImage.getBBox().getDimensions()
248  if nx != ny:
249  raise RuntimeError("Raw and calib maskedImages are trimmed differently in X and Y.")
250  if nx % 2 != 0:
251  raise RuntimeError("Calibration maskedImage is trimmed unevenly in X.")
252  if nx < 0:
253  raise RuntimeError("Calibration maskedImage is larger than raw data.")
254 
255  nEdge = nx//2
256  if nEdge > 0:
257  replacementMaskedImage = rawMaskedImage[nEdge:-nEdge, nEdge:-nEdge, afwImage.LOCAL]
258  SourceDetectionTask.setEdgeBits(
259  rawMaskedImage,
260  replacementMaskedImage.getBBox(),
261  rawMaskedImage.getMask().getPlaneBitMask("EDGE")
262  )
263  else:
264  replacementMaskedImage = rawMaskedImage
265 
266  return replacementMaskedImage
267 
268 

◆ updateVariance()

def lsst.ip.isr.isrFunctions.updateVariance (   maskedImage,
  gain,
  readNoise 
)
Set the variance plane based on the image plane.

Parameters
----------
maskedImage : `lsst.afw.image.MaskedImage`
    Image to process.  The variance plane is modified.
gain : scalar
    The amplifier gain in electrons/ADU.
readNoise : scalar
    The amplifier read nmoise in ADU/pixel.

Definition at line 342 of file isrFunctions.py.

342 def updateVariance(maskedImage, gain, readNoise):
343  """Set the variance plane based on the image plane.
344 
345  Parameters
346  ----------
347  maskedImage : `lsst.afw.image.MaskedImage`
348  Image to process. The variance plane is modified.
349  gain : scalar
350  The amplifier gain in electrons/ADU.
351  readNoise : scalar
352  The amplifier read nmoise in ADU/pixel.
353  """
354  var = maskedImage.getVariance()
355  var[:] = maskedImage.getImage()
356  var /= gain
357  var += readNoise**2
358 
359 

◆ widenSaturationTrails()

def lsst.ip.isr.isrFunctions.widenSaturationTrails (   mask)
Grow the saturation trails by an amount dependent on the width of the trail.

Parameters
----------
mask : `lsst.afw.image.Mask`
    Mask which will have the saturated areas grown.

Definition at line 775 of file isrFunctions.py.

775 def widenSaturationTrails(mask):
776  """Grow the saturation trails by an amount dependent on the width of the trail.
777 
778  Parameters
779  ----------
780  mask : `lsst.afw.image.Mask`
781  Mask which will have the saturated areas grown.
782  """
783 
784  extraGrowDict = {}
785  for i in range(1, 6):
786  extraGrowDict[i] = 0
787  for i in range(6, 8):
788  extraGrowDict[i] = 1
789  for i in range(8, 10):
790  extraGrowDict[i] = 3
791  extraGrowMax = 4
792 
793  if extraGrowMax <= 0:
794  return
795 
796  saturatedBit = mask.getPlaneBitMask("SAT")
797 
798  xmin, ymin = mask.getBBox().getMin()
799  width = mask.getWidth()
800 
801  thresh = afwDetection.Threshold(saturatedBit, afwDetection.Threshold.BITMASK)
802  fpList = afwDetection.FootprintSet(mask, thresh).getFootprints()
803 
804  for fp in fpList:
805  for s in fp.getSpans():
806  x0, x1 = s.getX0(), s.getX1()
807 
808  extraGrow = extraGrowDict.get(x1 - x0 + 1, extraGrowMax)
809  if extraGrow > 0:
810  y = s.getY() - ymin
811  x0 -= xmin + extraGrow
812  x1 -= xmin - extraGrow
813 
814  if x0 < 0:
815  x0 = 0
816  if x1 >= width - 1:
817  x1 = width - 1
818 
819  mask.array[y, x0:x1+1] |= saturatedBit
820 
821 
lsst::ip::isr.isrFunctions.interpolateDefectList
def interpolateDefectList(maskedImage, defectList, fwhm, fallbackValue=None)
Definition: isrFunctions.py:76
lsst::ip::isr.isrFunctions.applyGains
def applyGains(exposure, normalizeGains=False)
Definition: isrFunctions.py:745
lsst::ip::isr.isrFunctions.illuminationCorrection
def illuminationCorrection(maskedImage, illumMaskedImage, illumScale, trimToFit=True)
Definition: isrFunctions.py:410
lsst::ip::isr.isrFunctions.interpolateFromMask
def interpolateFromMask(maskedImage, fwhm, growSaturatedFootprints=1, maskNameList=['SAT'], fallbackValue=None)
Definition: isrFunctions.py:156
lsst::ip::isr.isrFunctions.saturationCorrection
def saturationCorrection(maskedImage, saturation, fwhm, growFootprints=1, interpolate=True, maskName='SAT', fallbackValue=None)
Definition: isrFunctions.py:189
lsst::ip::isr.isrFunctions.gainContext
def gainContext(exp, image, apply, gains=None)
Definition: isrFunctions.py:641
lsst::ip::isr.isrFunctions.transposeMaskedImage
def transposeMaskedImage(maskedImage)
Definition: isrFunctions.py:56
lsst::afw::math::FixedKernel
A kernel created from an Image.
Definition: Kernel.h:518
lsst::afw::math::makeStatistics
Statistics makeStatistics(lsst::afw::math::MaskedVector< EntryT > const &mv, std::vector< WeightPixel > const &vweights, int const flags, StatisticsControl const &sctrl=StatisticsControl())
The makeStatistics() overload to handle lsst::afw::math::MaskedVector<>
Definition: Statistics.h:520
lsst::afw::math::stringToStatisticsProperty
Property stringToStatisticsProperty(std::string const property)
Conversion function to switch a string to a Property (see Statistics.h)
Definition: Statistics.cc:747
lsst::ip::isr.isrFunctions.flatCorrection
def flatCorrection(maskedImage, flatMaskedImage, scalingType, userScale=1.0, invert=False, trimToFit=False)
Definition: isrFunctions.py:360
lsst::ip::isr.isrFunctions.darkCorrection
def darkCorrection(maskedImage, darkMaskedImage, expScale, darkScale, invert=False, trimToFit=False)
Definition: isrFunctions.py:298
lsst::ip::isr.isrFunctions.overscanCorrection
def overscanCorrection(ampMaskedImage, overscanImage, fitType='MEDIAN', order=1, collapseRej=3.0, statControl=None, overscanIsInt=True)
Definition: isrFunctions.py:441
lsst::ip::isr.isrFunctions.createPsf
def createPsf(fwhm)
Definition: isrFunctions.py:39
lsst::ip::isr.isrFunctions.brighterFatterCorrection
def brighterFatterCorrection(exposure, kernel, maxIter, threshold, applyGain, gains=None)
Definition: isrFunctions.py:518
lsst::ip::isr.isrFunctions.setBadRegions
def setBadRegions(exposure, badStatistic="MEDIAN")
Definition: isrFunctions.py:822
lsst::ip::isr.isrFunctions.growMasks
def growMasks(mask, radius=0, maskNameList=['BAD'], maskValue="BAD")
Definition: isrFunctions.py:135
lsst::ip::isr.isrFunctions.trimToMatchCalibBBox
def trimToMatchCalibBBox(rawMaskedImage, calibMaskedImage)
Definition: isrFunctions.py:222
lsst::ip::isr.isrFunctions.updateVariance
def updateVariance(maskedImage, gain, readNoise)
Definition: isrFunctions.py:342
lsst::afw::math::StatisticsControl
Pass parameters to a Statistics object.
Definition: Statistics.h:93
lsst::afw::math::ConvolutionControl
Parameters to control convolution.
Definition: ConvolveImage.h:50
lsst::afw::math::convolve
void convolve(OutImageT &convolvedImage, InImageT const &inImage, KernelT const &kernel, bool doNormalize, bool doCopyEdge=false)
Old, deprecated version of convolve.
Definition: ConvolveImage.cc:199
lsst::ip::isr.isrFunctions.attachTransmissionCurve
def attachTransmissionCurve(exposure, opticsTransmission=None, filterTransmission=None, sensorTransmission=None, atmosphereTransmission=None)
Definition: isrFunctions.py:693
lsst::ip::isr.isrFunctions.biasCorrection
def biasCorrection(maskedImage, biasMaskedImage, trimToFit=False)
Definition: isrFunctions.py:269
lsst::geom::Extent< int, 2 >
lsst::ip::isr.isrFunctions.makeThresholdMask
def makeThresholdMask(maskedImage, threshold, growFootprints=1, maskName='SAT')
Definition: isrFunctions.py:100
lsst::ip::isr.isrFunctions.widenSaturationTrails
def widenSaturationTrails(mask)
Definition: isrFunctions.py:775