24from scipy
import ndimage
29__all__ = [
"DcrModel",
"applyDcr",
"calculateDcr",
"calculateImageParallacticAngle"]
33 """A model of the true sky after correcting chromatic effects.
37 dcrNumSubfilters : `int`
38 Number of sub-filters used to model chromatic effects within a band.
40 A list of masked images, each containing the model for one subfilter
44 The ``DcrModel`` contains an estimate of the true sky, at a higher
45 wavelength resolution than the input observations. It can be forward-
46 modeled to produce Differential Chromatic Refraction (DCR) matched
47 templates
for a given ``Exposure``,
and provides utilities
for conditioning
48 the model
in ``dcrAssembleCoadd`` to avoid oscillating solutions between
49 iterations of forward modeling
or between the subfilters of the model.
52 def __init__(self, modelImages, effectiveWavelength, bandwidth, filterLabel=None, psf=None,
53 bbox=None, wcs=None, mask=None, variance=None, photoCalib=None):
67 def fromImage(cls, maskedImage, dcrNumSubfilters, effectiveWavelength, bandwidth,
68 wcs=None, filterLabel=None, psf=None, photoCalib=None):
69 """Initialize a DcrModel by dividing a coadd between the subfilters.
74 Input coadded image to divide equally between the subfilters.
75 dcrNumSubfilters : `int`
76 Number of sub-filters used to model chromatic effects within a
78 effectiveWavelength : `float`
79 The effective wavelengths of the current filter, in nanometers.
81 The bandwidth of the current filter,
in nanometers.
83 Coordinate system definition (wcs)
for the exposure.
85 The filter label, set
in the current instruments
' obs package.
86 Required for any calculation of DCR, including making matched
89 Point spread function (PSF) of the model.
90 Required
if the ``DcrModel`` will be persisted.
92 Calibration to convert instrumental flux
and
93 flux error to nanoJansky.
97 dcrModel : `lsst.pipe.tasks.DcrModel`
98 Best fit model of the true sky after correcting chromatic effects.
102 model = maskedImage.image.clone()
103 mask = maskedImage.mask.clone()
104 bbox = maskedImage.getBBox()
110 variance = maskedImage.variance.clone()
111 variance /= dcrNumSubfilters
112 model /= dcrNumSubfilters
113 modelImages = [model, ]
114 for subfilter
in range(1, dcrNumSubfilters):
115 modelImages.append(model.clone())
116 return cls(modelImages, effectiveWavelength, bandwidth,
117 filterLabel=filterLabel, psf=psf, bbox=bbox, wcs=wcs,
118 mask=mask, variance=variance, photoCalib=photoCalib)
121 def fromQuantum(cls, availableCoaddRefs, effectiveWavelength, bandwidth, numSubfilters):
122 """Load an existing DcrModel from a Gen 3 repository.
126 availableCoaddRefs : `dict` [`int`, `lsst.daf.butler.DeferredDatasetHandle`]
127 Dictionary of spatially relevant retrieved coadd patches,
128 indexed by their sequential patch number.
129 effectiveWavelength : `float`
130 The effective wavelengths of the current filter, in nanometers.
132 The bandwidth of the current filter,
in nanometers.
133 numSubfilters : `int`
134 Number of subfilters
in the DcrCoadd.
138 dcrModel : `lsst.pipe.tasks.DcrModel`
139 Best fit model of the true sky after correcting chromatic effects.
148 modelImages = [
None]*numSubfilters
150 for coaddRef
in availableCoaddRefs:
151 subfilter = coaddRef.dataId[
"subfilter"]
152 dcrCoadd = coaddRef.get()
153 if filterLabel
is None:
154 filterLabel = dcrCoadd.getFilter()
156 psf = dcrCoadd.getPsf()
158 bbox = dcrCoadd.getBBox()
164 variance = dcrCoadd.variance
165 if photoCalib
is None:
166 photoCalib = dcrCoadd.getPhotoCalib()
167 modelImages[subfilter] = dcrCoadd.image
168 return cls(modelImages, effectiveWavelength, bandwidth, filterLabel,
169 psf, bbox, wcs, mask, variance, photoCalib)
172 """Return the number of subfilters.
176 dcrNumSubfilters : `int`
177 The number of DCR subfilters in the model.
182 """Iterate over the subfilters of the DCR model.
187 Index of the current ``subfilter`` within the full band.
188 Negative indices are allowed, and count
in reverse order
189 from the highest ``subfilter``.
194 The DCR model
for the given ``subfilter``.
199 If the requested ``subfilter``
is greater
or equal to the number
200 of subfilters
in the model.
202 if np.abs(subfilter) >= len(self):
203 raise IndexError(
"subfilter out of bounds.")
207 """Update the model image for one subfilter.
212 Index of the current subfilter within the full band.
214 The DCR model to set for the given ``subfilter``.
219 If the requested ``subfilter``
is greater
or equal to the number
220 of subfilters
in the model.
222 If the bounding box of the new image does
not match.
224 if np.abs(subfilter) >= len(self):
225 raise IndexError(
"subfilter out of bounds.")
226 if maskedImage.getBBox() != self.
bbox:
227 raise ValueError(
"The bounding box of a subfilter must not change.")
232 """Return the effective wavelength of the model.
236 effectiveWavelength : `float`
237 The effective wavelength of the current filter, in nanometers.
243 """Return the filter label for the model.
248 The filter used for the input observations.
254 """Return the bandwidth of the model.
259 The bandwidth of the current filter, in nanometers.
265 """Return the psf of the model.
270 Point spread function (PSF) of the model.
276 """Return the common bounding box of each subfilter image.
280 bbox : `lsst.afw.geom.Box2I`
281 Bounding box of the DCR model.
287 """Return the WCS of each subfilter image.
292 Coordinate system definition (wcs) for the exposure.
298 """Return the common mask of each subfilter image.
303 Mask plane of the DCR model.
309 """Return the common variance of each subfilter image.
314 Variance plane of the DCR model.
319 """Calculate a reference image from the average of the subfilter
324 bbox : `lsst.afw.geom.Box2I`, optional
325 Sub-region of the coadd. Returns the entire image if `
None`.
329 refImage : `numpy.ndarray`
330 The reference image
with no chromatic effects applied.
332 bbox = bbox or self.
bbox
333 return np.mean([model[bbox].array
for model
in self], axis=0)
335 def assign(self, dcrSubModel, bbox=None):
336 """Update a sub-region of the ``DcrModel`` with new values.
340 dcrSubModel : `lsst.pipe.tasks.DcrModel`
341 New model of the true scene after correcting chromatic effects.
342 bbox : `lsst.afw.geom.Box2I`, optional
343 Sub-region of the coadd.
344 Defaults to the bounding box of ``dcrSubModel``.
349 If the new model has a different number of subfilters.
351 if len(dcrSubModel) != len(self):
352 raise ValueError(
"The number of DCR subfilters must be the same "
353 "between the old and new models.")
354 bbox = bbox
or self.
bbox
355 for model, subModel
in zip(self, dcrSubModel):
356 model.assign(subModel[bbox], bbox)
359 visitInfo=None, bbox=None, mask=None,
360 splitSubfilters=True, splitThreshold=0., amplifyModel=1.):
361 """Create a DCR-matched template image for an exposure.
366 The input exposure to build a matched template for.
367 May be omitted
if all of the metadata
is supplied separately
368 order : `int`, optional
369 Interpolation order of the DCR shift.
371 Metadata
for the exposure. Ignored
if ``exposure``
is set.
372 bbox : `lsst.afw.geom.Box2I`, optional
373 Sub-region of the coadd,
or use the entire coadd
if not supplied.
375 reference mask to use
for the template image.
376 splitSubfilters : `bool`, optional
377 Calculate DCR
for two evenly-spaced wavelengths
in each subfilter,
378 instead of at the midpoint. Default:
True
379 splitThreshold : `float`, optional
380 Minimum DCR difference within a subfilter required to use
382 amplifyModel : `float`, optional
383 Multiplication factor to amplify differences between model planes.
384 Used to speed convergence of iterative forward modeling.
388 templateImage : `lsst.afw.image.ImageF`
389 The DCR-matched template
394 If neither ``exposure``
or ``visitInfo`` are set.
397 raise ValueError(
"'effectiveWavelength' and 'bandwidth' must be set for the DcrModel in order "
399 if exposure
is not None:
400 visitInfo = exposure.getInfo().getVisitInfo()
401 elif visitInfo
is None:
402 raise ValueError(
"Either exposure or visitInfo must be set.")
406 splitSubfilters=splitSubfilters)
407 templateImage = afwImage.ImageF(bbox)
409 for subfilter, dcr
in enumerate(dcrShift):
410 if self[subfilter]
is None:
412 self.log.debug(
"Skipping missing DCR model subfilter %d", subfilter)
419 model = (self[subfilter][bbox].array - refModel)*amplifyModel + refModel
421 model = self[subfilter][bbox].array
422 templateImage.array += applyDcr(model, dcr, splitSubfilters=splitSubfilters,
423 splitThreshold=splitThreshold, order=order)
427 visitInfo=None, bbox=None, mask=None):
428 """Wrapper to create an exposure from a template image.
433 The input exposure to build a matched template for.
434 May be omitted
if all of the metadata
is supplied separately
436 Metadata
for the exposure. Ignored
if ``exposure``
is set.
437 bbox : `lsst.afw.geom.Box2I`, optional
438 Sub-region of the coadd,
or use the entire coadd
if not supplied.
440 reference mask to use
for the template image.
444 templateExposure : `lsst.afw.image.exposureF`
445 The DCR-matched template
450 If no `photcCalib`
is set.
455 bbox=bbox, mask=mask)
456 maskedImage = afwImage.MaskedImageF(bbox)
457 maskedImage.image = templateImage[bbox]
458 maskedImage.mask = self.
mask[bbox]
459 maskedImage.variance = self.
variance[bbox]
463 templateExposure = afwImage.ExposureF(bbox, self.
wcs)
464 templateExposure.setMaskedImage(maskedImage[bbox])
465 templateExposure.setPsf(self.
psf)
466 templateExposure.setFilter(self.
filter)
468 raise RuntimeError(
"No PhotoCalib set for the DcrModel. "
469 "If the DcrModel was created from a masked image"
470 " you must also specify the photoCalib.")
471 templateExposure.setPhotoCalib(self.
photoCalib)
472 return templateExposure
475 """Average two iterations' solutions to reduce oscillations.
480 The new DCR model images from the current iteration.
481 The values will be modified
in place.
482 bbox : `lsst.afw.geom.Box2I`
483 Sub-region of the coadd
484 gain : `float`, optional
485 Relative weight to give the new solution when updating the model.
486 Defaults to 1.0, which gives equal weight to both solutions.
489 for model, newModel
in zip(self, modelImages):
491 newModel += model[bbox]
492 newModel /= 1. + gain
495 regularizationWidth=2):
496 """Restrict large variations in the model between iterations.
501 Index of the current subfilter within the full band.
503 The new DCR model for one subfilter
from the current iteration.
504 Values
in ``newModel`` that are extreme compared
with the last
505 iteration are modified
in place.
506 bbox : `lsst.afw.geom.Box2I`
508 regularizationFactor : `float`
509 Maximum relative change of the model allowed between iterations.
510 regularizationWidth : int, optional
511 Minimum radius of a region to include
in regularization,
in pixels.
513 refImage = self[subfilter][bbox].array
514 highThreshold = np.abs(refImage)*regularizationFactor
515 lowThreshold = refImage/regularizationFactor
516 newImage = newModel.array
518 regularizationWidth=regularizationWidth)
521 regularizationWidth=2, mask=None, convergenceMaskPlanes="DETECTED"):
522 """Restrict large variations in the model between subfilters.
527 The new DCR model images from the current iteration.
528 The values will be modified
in place.
529 bbox : `lsst.afw.geom.Box2I`
532 Statistics control object
for coaddition.
533 regularizationFactor : `float`
534 Maximum relative change of the model allowed between subfilters.
535 regularizationWidth : `int`, optional
536 Minimum radius of a region to include
in regularization,
in pixels.
538 Optional alternate mask
539 convergenceMaskPlanes : `list` of `str`,
or `str`, optional
540 Mask planes to use to calculate convergence.
544 This implementation of frequency regularization restricts each
545 subfilter image to be a smoothly-varying function times a reference
551 maxDiff = np.sqrt(regularizationFactor)
552 noiseLevel = self.
calculateNoiseCutoff(modelImages[0], statsCtrl, bufferSize=5, mask=mask, bbox=bbox)
554 badPixels = np.isnan(referenceImage) | (referenceImage <= 0.)
555 if np.sum(~badPixels) == 0:
558 referenceImage[badPixels] = 0.
559 filterWidth = regularizationWidth
560 fwhm = 2.*filterWidth
564 smoothRef = ndimage.gaussian_filter(referenceImage, filterWidth, mode=
'constant')
568 smoothRef += 3.*noiseLevel
570 lowThreshold = smoothRef/maxDiff
571 highThreshold = smoothRef*maxDiff
572 for model
in modelImages:
574 highThreshold=highThreshold,
575 lowThreshold=lowThreshold,
576 regularizationWidth=regularizationWidth)
577 smoothModel = ndimage.gaussian_filter(model.array, filterWidth, mode=
'constant')
578 smoothModel += 3.*noiseLevel
579 relativeModel = smoothModel/smoothRef
582 relativeModel2 = ndimage.gaussian_filter(relativeModel, filterWidth/alpha)
583 relativeModel += alpha*(relativeModel - relativeModel2)
584 model.array = relativeModel*referenceImage
587 convergenceMaskPlanes="DETECTED", mask=None, bbox=None):
588 """Helper function to calculate the background noise level of an image.
593 The input image to evaluate the background noise properties.
595 Statistics control object for coaddition.
597 Number of additional pixels to exclude
598 from the edges of the bounding box.
599 convergenceMaskPlanes : `list` of `str`,
or `str`
600 Mask planes to use to calculate convergence.
602 Optional alternate mask
603 bbox : `lsst.afw.geom.Box2I`, optional
604 Sub-region of the masked image to calculate the noise level over.
608 noiseCutoff : `float`
609 The threshold value to treat pixels
as noise
in an image..
614 mask = self.
mask[bbox]
616 bboxShrink.grow(-bufferSize)
617 convergeMask = mask.getPlaneBitMask(convergenceMaskPlanes)
619 backgroundPixels = mask[bboxShrink].array & (statsCtrl.getAndMask() | convergeMask) == 0
620 noiseCutoff = np.std(image[bboxShrink].array[backgroundPixels])
624 """Restrict image values to be between upper and lower limits.
626 This method flags all pixels in an image that are outside of the given
627 threshold values. The threshold values are taken
from a reference
628 image, so noisy pixels are likely to get flagged. In order to exclude
629 those noisy pixels, the array of flags
is eroded
and dilated, which
630 removes isolated pixels outside of the thresholds
from the list of
631 pixels to be modified. Pixels that remain flagged after this operation
632 have their values set to the appropriate upper
or lower threshold
637 image : `numpy.ndarray`
638 The image to apply the thresholds to.
639 The values will be modified
in place.
640 highThreshold : `numpy.ndarray`, optional
641 Array of upper limit values
for each pixel of ``image``.
642 lowThreshold : `numpy.ndarray`, optional
643 Array of lower limit values
for each pixel of ``image``.
644 regularizationWidth : `int`, optional
645 Minimum radius of a region to include
in regularization,
in pixels.
650 filterStructure = ndimage.iterate_structure(ndimage.generate_binary_structure(2, 1),
652 if highThreshold
is not None:
653 highPixels = image > highThreshold
654 if regularizationWidth > 0:
656 highPixels = ndimage.binary_opening(highPixels, structure=filterStructure)
657 image[highPixels] = highThreshold[highPixels]
658 if lowThreshold
is not None:
659 lowPixels = image < lowThreshold
660 if regularizationWidth > 0:
662 lowPixels = ndimage.binary_opening(lowPixels, structure=filterStructure)
663 image[lowPixels] = lowThreshold[lowPixels]
666def applyDcr(image, dcr, useInverse=False, splitSubfilters=False, splitThreshold=0.,
667 doPrefilter=True, order=3):
668 """Shift an image along the X and Y directions.
672 image : `numpy.ndarray`
673 The input image to shift.
675 Shift calculated with ``calculateDcr``.
676 Uses numpy axes ordering (Y, X).
677 If ``splitSubfilters``
is set, each element
is itself a `tuple`
678 of two `float`, corresponding to the DCR shift at the two wavelengths.
679 Otherwise, each element
is a `float` corresponding to the DCR shift at
680 the effective wavelength of the subfilter.
681 useInverse : `bool`, optional
682 Apply the shift
in the opposite direction. Default:
False
683 splitSubfilters : `bool`, optional
684 Calculate DCR
for two evenly-spaced wavelengths
in each subfilter,
685 instead of at the midpoint. Default:
False
686 splitThreshold : `float`, optional
687 Minimum DCR difference within a subfilter required to use
689 doPrefilter : `bool`, optional
690 Spline filter the image before shifting,
if set. Filtering
is required,
691 so only set to
False if the image
is already filtered.
692 Filtering takes ~20% of the time of shifting, so
if `applyDcr` will be
693 called repeatedly on the same image it
is more efficient to
694 precalculate the filter.
695 order : `int`, optional
696 The order of the spline interpolation, default
is 3.
700 shiftedImage : `numpy.ndarray`
701 A copy of the input image
with the specified shift applied.
703 if doPrefilter
and order > 1:
704 prefilteredImage = ndimage.spline_filter(image, order=order)
706 prefilteredImage = image
708 shiftAmp = np.max(np.abs([_dcr0 - _dcr1
for _dcr0, _dcr1
in zip(dcr[0], dcr[1])]))
709 if shiftAmp >= splitThreshold:
711 shift = [-1.*s
for s
in dcr[0]]
712 shift1 = [-1.*s
for s
in dcr[1]]
716 shiftedImage = ndimage.shift(prefilteredImage, shift, prefilter=
False, order=order)
717 shiftedImage += ndimage.shift(prefilteredImage, shift1, prefilter=
False, order=order)
723 dcr = (np.mean(dcr[0]), np.mean(dcr[1]))
725 shift = [-1.*s
for s
in dcr]
728 shiftedImage = ndimage.shift(prefilteredImage, shift, prefilter=
False, order=order)
732def calculateDcr(visitInfo, wcs, effectiveWavelength, bandwidth, dcrNumSubfilters, splitSubfilters=False):
733 """Calculate the shift in pixels of an exposure due to DCR.
738 Metadata for the exposure.
740 Coordinate system definition (wcs)
for the exposure.
741 effectiveWavelength : `float`
742 The effective wavelengths of the current filter,
in nanometers.
744 The bandwidth of the current filter,
in nanometers.
745 dcrNumSubfilters : `int`
746 Number of sub-filters used to model chromatic effects within a band.
747 splitSubfilters : `bool`, optional
748 Calculate DCR
for two evenly-spaced wavelengths
in each subfilter,
749 instead of at the midpoint. Default:
False
753 dcrShift : `tuple` of two `float`
754 The 2D shift due to DCR,
in pixels.
755 Uses numpy axes ordering (Y, X).
759 weight = [0.75, 0.25]
763 diffRefractAmp0 = differentialRefraction(wavelength=wl0, wavelengthRef=effectiveWavelength,
764 elevation=visitInfo.getBoresightAzAlt().getLatitude(),
765 observatory=visitInfo.getObservatory(),
766 weather=visitInfo.getWeather())
767 diffRefractAmp1 = differentialRefraction(wavelength=wl1, wavelengthRef=effectiveWavelength,
768 elevation=visitInfo.getBoresightAzAlt().getLatitude(),
769 observatory=visitInfo.getObservatory(),
770 weather=visitInfo.getWeather())
772 diffRefractPix0 = diffRefractAmp0.asArcseconds()/wcs.getPixelScale().asArcseconds()
773 diffRefractPix1 = diffRefractAmp1.asArcseconds()/wcs.getPixelScale().asArcseconds()
774 diffRefractArr = [diffRefractPix0*weight[0] + diffRefractPix1*weight[1],
775 diffRefractPix0*weight[1] + diffRefractPix1*weight[0]]
776 shiftX = [diffRefractPix*np.sin(rotation.asRadians())
for diffRefractPix
in diffRefractArr]
777 shiftY = [diffRefractPix*np.cos(rotation.asRadians())
for diffRefractPix
in diffRefractArr]
778 dcrShift.append(((shiftY[0], shiftX[0]), (shiftY[1], shiftX[1])))
780 diffRefractAmp = (diffRefractAmp0 + diffRefractAmp1)/2.
781 diffRefractPix = diffRefractAmp.asArcseconds()/wcs.getPixelScale().asArcseconds()
782 shiftX = diffRefractPix*np.sin(rotation.asRadians())
783 shiftY = diffRefractPix*np.cos(rotation.asRadians())
784 dcrShift.append((shiftY, shiftX))
789 """Calculate the total sky rotation angle of an exposure.
794 Metadata for the exposure.
796 Coordinate system definition (wcs)
for the exposure.
801 The rotation of the image axis, East
from North.
802 Equal to the parallactic angle plus any additional rotation of the
804 A rotation angle of 0 degrees
is defined
with
805 North along the +y axis
and East along the +x axis.
806 A rotation angle of 90 degrees
is defined
with
807 North along the +x axis
and East along the -y axis.
809 parAngle = visitInfo.getBoresightParAngle().asRadians()
810 cd = wcs.getCdMatrix()
812 cdAngle = (np.arctan2(-cd[0, 1], cd[0, 0]) + np.arctan2(cd[1, 0], cd[1, 1]))/2.
813 rotAngle = (cdAngle + parAngle)*geom.radians
815 cdAngle = (np.arctan2(cd[0, 1], -cd[0, 0]) + np.arctan2(cd[1, 0], cd[1, 1]))/2.
816 rotAngle = (cdAngle - parAngle)*geom.radians
821 """Iterate over the wavelength endpoints of subfilters.
825 effectiveWavelength : `float`
826 The effective wavelength of the current filter, in nanometers.
828 The bandwidth of the current filter,
in nanometers.
829 dcrNumSubfilters : `int`
830 Number of sub-filters used to model chromatic effects within a band.
834 `tuple` of two `float`
835 The next set of wavelength endpoints
for a subfilter,
in nanometers.
837 lambdaMin = effectiveWavelength - bandwidth/2
838 lambdaMax = effectiveWavelength + bandwidth/2
839 wlStep = bandwidth/dcrNumSubfilters
840 for wl
in np.linspace(lambdaMin, lambdaMax, dcrNumSubfilters, endpoint=
False):
841 yield (wl, wl + wlStep)
table::Key< table::Array< std::uint8_t > > wcs
A polymorphic base class for representing an image's Point Spread Function.
A 2-dimensional celestial WCS that transform pixels to ICRS RA/Dec, using the LSST standard for pixel...
A class to contain the data, WCS, and other information needed to describe an image of the sky.
A group of labels for a filter in an exposure or coadd.
A class to represent a 2-dimensional array of pixels.
Represent a 2-dimensional array of bitmask pixels.
A class to manipulate images, masks, and variance as a single object.
The photometric calibration of an exposure.
Information about a single exposure of an imaging camera.
Pass parameters to a Statistics object.
A class representing an angle.
An integer coordinate rectangle.
def fromImage(cls, maskedImage, dcrNumSubfilters, effectiveWavelength, bandwidth, wcs=None, filterLabel=None, psf=None, photoCalib=None)
def calculateNoiseCutoff(self, image, statsCtrl, bufferSize, convergenceMaskPlanes="DETECTED", mask=None, bbox=None)
def applyImageThresholds(self, image, highThreshold=None, lowThreshold=None, regularizationWidth=2)
def __init__(self, modelImages, effectiveWavelength, bandwidth, filterLabel=None, psf=None, bbox=None, wcs=None, mask=None, variance=None, photoCalib=None)
def getReferenceImage(self, bbox=None)
def conditionDcrModel(self, modelImages, bbox, gain=1.)
def regularizeModelIter(self, subfilter, newModel, bbox, regularizationFactor, regularizationWidth=2)
def assign(self, dcrSubModel, bbox=None)
def regularizeModelFreq(self, modelImages, bbox, statsCtrl, regularizationFactor, regularizationWidth=2, mask=None, convergenceMaskPlanes="DETECTED")
def __setitem__(self, subfilter, maskedImage)
def buildMatchedExposure(self, exposure=None, visitInfo=None, bbox=None, mask=None)
def __getitem__(self, subfilter)
def fromQuantum(cls, availableCoaddRefs, effectiveWavelength, bandwidth, numSubfilters)
def buildMatchedTemplate(self, exposure=None, order=3, visitInfo=None, bbox=None, mask=None, splitSubfilters=True, splitThreshold=0., amplifyModel=1.)
def effectiveWavelength(self)
def wavelengthGenerator(effectiveWavelength, bandwidth, dcrNumSubfilters)
def calculateImageParallacticAngle(visitInfo, wcs)