23 from contextlib
import contextmanager
31 from .anetBasicAstrometry
import ANetBasicAstrometryTask
32 from .sip
import makeCreateWcsWithSip
35 solver = pexConfig.ConfigField(
36 dtype = ANetBasicAstrometryTask.ConfigClass,
37 doc =
"Configuration for the astrometry solver",
39 forceKnownWcs = pexConfig.Field(dtype=bool, doc=(
40 "Assume that the input image's WCS is correct, without comparing it to any external reality." +
41 " (In contrast to using Astrometry.net). NOTE, if you set this, you probably also want to" +
42 " un-set 'solver.calculateSip'; otherwise we'll still try to find a TAN-SIP WCS starting " +
43 " from the existing WCS"), default=
False)
44 rejectThresh = pexConfig.RangeField(dtype=float, default=3.0, doc=
"Rejection threshold for Wcs fitting",
45 min=0.0, inclusiveMin=
False)
46 rejectIter = pexConfig.RangeField(dtype=int, default=3, doc=
"Rejection iterations for Wcs fitting",
58 """!Use astrometry.net to match input sources with a reference catalog and solve for the Wcs
60 The actual matching and solving is done by the 'solver'; this Task
61 serves as a wrapper for taking into account the known optical distortion.
63 \section pipe_tasks_astrometry_Contents Contents
65 - \ref pipe_tasks_astrometry_Purpose
66 - \ref pipe_tasks_astrometry_Initialize
67 - \ref pipe_tasks_astrometry_IO
68 - \ref pipe_tasks_astrometry_Config
69 - \ref pipe_tasks_astrometry_Debug
70 - \ref pipe_tasks_astrometry_Example
72 \section pipe_tasks_astrometry_Purpose Description
74 \copybrief ANetAstrometryTask
76 \section pipe_tasks_astrometry_Initialize Task initialisation
80 \section pipe_tasks_astrometry_IO Invoking the Task
84 \section pipe_tasks_astrometry_Config Configuration parameters
86 See \ref ANetAstrometryConfig
88 \section pipe_tasks_astrometry_Debug Debug variables
90 The \link lsst.pipe.base.cmdLineTask.CmdLineTask command line task\endlink interface supports a
91 flag \c -d to import \b debug.py from your \c PYTHONPATH; see \ref baseDebug for more about \b debug.py files.
93 The available variables in ANetAstrometryTask are:
96 <DD> If True call astrometry.showAstrometry while iterating ANetAstrometryConfig.rejectIter times,
97 and also after converging.
99 <DD> ds9 frame to use in astrometry.showAstrometry
101 <DD> Pause within astrometry.showAstrometry
104 \section pipe_tasks_astrometry_Example A complete example of using ANetAstrometryTask
106 See \ref meas_photocal_photocal_Example.
108 To investigate the \ref pipe_tasks_astrometry_Debug, put something like
112 di = lsstDebug.getInfo(name) # N.b. lsstDebug.Info(name) would call us recursively
113 if name == "lsst.pipe.tasks.astrometry":
118 lsstDebug.Info = DebugInfo
120 into your debug.py file and run photoCalTask.py with the \c --debug flag.
122 ConfigClass = ANetAstrometryConfig
125 """!Create the astrometric calibration task. Most arguments are simply passed onto pipe.base.Task.
127 \param schema An lsst::afw::table::Schema used to create the output lsst.afw.table.SourceCatalog
128 \param **kwds keyword arguments to be passed to the lsst.pipe.base.task.Task constructor
130 A centroid field "centroid.distorted" (used internally during the Task's operation)
131 will be added to the schema.
133 pipeBase.Task.__init__(self, **kwds)
135 if schema.getVersion() == 0:
138 doc=
"centroid distorted for astrometry solver")
140 doc=
"centroid distorted err for astrometry solver")
142 doc=
"centroid distorted flag astrometry solver")
146 doc=
"centroid distorted for astrometry solver")
148 doc=
"centroid distorted for astrometry solver")
150 doc=
"centroid distorted err for astrometry solver")
152 doc=
"centroid distorted err for astrometry solver")
154 doc=
"centroid distorted flag astrometry solver")
161 def run(self, exposure, sourceCat):
162 """!Match with reference sources and calculate an astrometric solution
164 \param[in,out] exposure Exposure to calibrate
165 \param[in] sourceCat catalog of measured sources (an lsst.afw.table.SourceCatalog)
166 \return a pipeBase.Struct with fields:
167 - matches: Astrometric matches, as an lsst.afw.table.ReferenceMatchVector
168 - matchMeta: Metadata for astrometric matches
170 The reference catalog actually used is up to the implementation
171 of the solver; it will be manifested in the returned matches as
172 a list of lsst.afw.table.ReferenceMatch objects (\em i.e. of lsst.afw.table.Match with
173 \c first being of type lsst.afw.table.SimpleRecord and \c second type lsst.afw.table.SourceRecord ---
174 the reference object and matched object respectively).
177 The input sources have the centroid slot moved to a new column "centroid.distorted"
178 which has the positions corrected for any known optical distortion;
179 the 'solver' (which is instantiated in the 'astrometry' member)
180 should therefore simply use the centroids provided by calling
181 afw.table.Source.getCentroid() on the individual source records. This column \em must
182 be present in the sources table.
185 results = self.
astrometry(sourceCat=sourceCat, exposure=exposure, bbox=bbox)
188 self.
refitWcs(sourceCat=sourceCat, exposure=exposure, matches=results.matches)
194 """!Calculate distorted source positions
196 CCD images are often affected by optical distortion that makes
197 the astrometric solution higher order than linear. Unfortunately,
198 most (all?) matching algorithms require that the distortion be
199 small or zero, and so it must be removed. We do this by calculating
200 (un-)distorted positions, based on a known optical distortion model
203 The distortion correction moves sources, so we return the distorted bounding box.
205 \param[in] exposure Exposure to process
206 \param[in,out] sourceCat SourceCatalog; getX() and getY() will be used as inputs,
207 with distorted points in "centroid.distorted" field.
208 \return bounding box of distorted exposure
210 detector = exposure.getDetector()
211 pixToTanXYTransform =
None
213 self.log.warn(
"No detector associated with exposure; assuming null distortion")
215 tanSys = detector.makeCameraSys(TAN_PIXELS)
216 pixToTanXYTransform = detector.getTransformMap().get(tanSys)
218 if pixToTanXYTransform
is None:
219 self.log.info(
"Null distortion correction")
224 return exposure.getBBox()
227 self.log.info(
"Applying distortion correction")
229 centroid = pixToTanXYTransform.forwardTransform(s.getCentroid())
237 for corner
in detector.getCorners(TAN_PIXELS):
238 bboxD.include(corner)
243 """!Context manager that applies and removes distortion
245 We move the "centroid" definition in the catalog table to
246 point to the distorted positions. This is undone on exit
249 The input Wcs is taken to refer to the coordinate system
250 with the distortion correction applied, and hence no shift
251 is required when the sources are distorted. However, after
252 Wcs fitting, the Wcs is in the distorted frame so when the
253 distortion correction is removed, the Wcs needs to be
254 shifted to compensate.
256 \param sourceCat Sources on exposure, an lsst.afw.table.SourceCatalog
257 \param exposure Exposure holding Wcs, an lsst.afw.image.ExposureF or D
258 \return bounding box of distorted exposure
261 bbox = self.
distort(sourceCat=sourceCat, exposure=exposure)
262 oldCentroidName = sourceCat.table.getCentroidDefinition()
268 sourceCat.table.defineCentroid(oldCentroidName)
269 x0, y0 = exposure.getXY0()
270 wcs = exposure.getWcs()
272 wcs.shiftReferencePixel(-bbox.getMinX() + x0, -bbox.getMinY() + y0)
276 """!Solve astrometry to produce WCS
278 \param sourceCat Sources on exposure, an lsst.afw.table.SourceCatalog
279 \param exposure Exposure to process, an lsst.afw.image.ExposureF or D
280 \param bbox Bounding box, or None to use exposure
281 \return a pipe.base.Struct with fields:
282 - matches: star matches
283 - matchMeta: match metadata
285 if not self.config.forceKnownWcs:
286 self.log.info(
"Solving astrometry")
288 bbox = exposure.getBBox()
291 self.
astrometer = ANetBasicAstrometryTask(self.config.solver, log=self.log)
293 if self.config.forceKnownWcs:
294 self.log.info(
"Forcing the input exposure's WCS")
295 if self.config.solver.calculateSip:
296 self.log.warn(
"'forceKnownWcs' and 'solver.calculateSip' options are both set." +
297 " Will try to compute a TAN-SIP WCS starting from the input WCS.")
298 astrom = self.astrometer.useKnownWcs(sourceCat=sourceCat, exposure=exposure, bbox=bbox)
300 astrom = self.astrometer.determineWcs(sourceCat=sourceCat, exposure=exposure, bbox=bbox)
302 if astrom
is None or astrom.getWcs()
is None:
303 raise RuntimeError(
"Unable to solve astrometry")
305 matches = astrom.getMatches()
306 matchMeta = astrom.getMatchMetadata()
307 if matches
is None or len(matches) == 0:
308 raise RuntimeError(
"No astrometric matches")
309 self.log.info(
"%d astrometric matches" % (len(matches)))
311 if not self.config.forceKnownWcs:
313 exposure.setWcs(astrom.getWcs())
315 self.display(
'astrometry', exposure=exposure, sources=sourceCat, matches=matches)
317 return pipeBase.Struct(matches=matches, matchMeta=matchMeta)
321 """!A final Wcs solution after matching and removing distortion
323 Specifically, fitting the non-linear part, since the linear
324 part has been provided by the matching engine.
326 \param sourceCat Sources on exposure, an lsst.afw.table.SourceCatalog
327 \param exposure Exposure of interest, an lsst.afw.image.ExposureF or D
328 \param matches Astrometric matches, as an lsst.afw.table.ReferenceMatchVector
330 \return the resolved-Wcs object, or None if config.solver.calculateSip is False.
333 if self.config.solver.calculateSip:
334 self.log.info(
"Refitting WCS")
335 origMatches = matches
336 wcs = exposure.getWcs()
343 def fitWcs(initialWcs, title=None):
344 """!Do the WCS fitting and display of the results"""
346 resultWcs = sip.getNewWcs()
348 showAstrometry(exposure, resultWcs, origMatches, matches, frame=frame,
349 title=title, pause=pause)
350 return resultWcs, sip.getScatterOnSky()
354 for i
in range(self.config.rejectIter):
355 wcs, scatter = fitWcs(wcs, title=
"Iteration %d" % i)
357 ref = numpy.array([wcs.skyToPixel(m.first.getCoord())
for m
in matches])
358 src = numpy.array([m.second.getCentroid()
for m
in matches])
362 for d, m
in zip(diff, matches):
363 if numpy.all(numpy.abs(d) < self.config.rejectThresh*rms):
367 if len(matches) == len(trimmed):
372 wcs, scatter = fitWcs(wcs, title=
"Final astrometry")
374 except lsst.pex.exceptions.LengthError
as e:
375 self.log.warn(
"Unable to fit SIP: %s" % e)
377 self.log.info(
"Astrometric scatter: %f arcsec (%s non-linear terms, %d matches, %d rejected)" %
378 (scatter.asArcseconds(),
"with" if wcs.hasDistortion()
else "without",
379 len(matches), numRejected))
383 for index, source
in enumerate(sourceCat):
384 sky = wcs.pixelToSky(source.getX(), source.getY())
387 self.log.warn(
"Not calculating a SIP solution; matches may be suspect")
389 self.display(
'astrometry', exposure=exposure, sources=sourceCat, matches=matches)
394 def showAstrometry(exposure, wcs, allMatches, useMatches, frame=0, title=None, pause=False):
395 """!Show results of astrometry fitting
397 \param exposure Image to display
398 \param wcs Astrometric solution
399 \param allMatches List of all astrometric matches (including rejects)
400 \param useMatches List of used astrometric matches
401 \param frame Frame number for display
402 \param title Title for display
403 \param pause Pause to allow viewing of the display and optional debugging?
405 - Matches are shown in yellow if used in the Wcs solution, otherwise red
406 - +: Detected objects
407 - x: Catalogue objects
410 ds9.mtv(exposure, frame=frame, title=title)
412 useIndices = set(m.second.getId()
for m
in useMatches)
415 with ds9.Buffering():
416 for i, m
in enumerate(allMatches):
417 x, y = m.second.getX(), m.second.getY()
418 pix = wcs.skyToPixel(m.first.getCoord())
420 isUsed = m.second.getId()
in useIndices
422 radii.append(numpy.hypot(pix[0] - x, pix[1] - y))
424 color = ds9.YELLOW
if isUsed
else ds9.RED
426 ds9.dot(
"+", x, y, size=10, frame=frame, ctype=color)
427 ds9.dot(
"x", pix[0], pix[1], size=10, frame=frame, ctype=color)
429 radii = numpy.array(radii)
430 print "<dr> = %.4g +- %.4g pixels [%d/%d matches]" % (radii.mean(), radii.std(),
431 len(useMatches), len(allMatches))
437 reply = raw_input(
"Debugging? [p]db [q]uit; any other key to continue... ").strip()
441 reply = reply.split()
445 import pdb;pdb.set_trace()
def refitWcs
A final Wcs solution after matching and removing distortion.
def __init__
Create the astrometric calibration task.
An integer coordinate rectangle.
CreateWcsWithSip< MatchT > makeCreateWcsWithSip(std::vector< MatchT > const &matches, afw::image::Wcs const &linearWcs, int const order, afw::geom::Box2I const &bbox=afw::geom::Box2I(), int const ngrid=0)
Factory function for CreateWcsWithSip.
def distortionContext
Context manager that applies and removes distortion.
Use astrometry.net to match input sources with a reference catalog and solve for the Wcs...
PointKey< double > Point2DKey
def astrometry
Solve astrometry to produce WCS.
A floating-point coordinate rectangle geometry.
def run
Match with reference sources and calculate an astrometric solution.
def distort
Calculate distorted source positions.
def showAstrometry
Show results of astrometry fitting.