LSST Applications 27.0.0,g0265f82a02+469cd937ee,g02d81e74bb+21ad69e7e1,g1470d8bcf6+cbe83ee85a,g2079a07aa2+e67c6346a6,g212a7c68fe+04a9158687,g2305ad1205+94392ce272,g295015adf3+81dd352a9d,g2bbee38e9b+469cd937ee,g337abbeb29+469cd937ee,g3939d97d7f+72a9f7b576,g487adcacf7+71499e7cba,g50ff169b8f+5929b3527e,g52b1c1532d+a6fc98d2e7,g591dd9f2cf+df404f777f,g5a732f18d5+be83d3ecdb,g64a986408d+21ad69e7e1,g858d7b2824+21ad69e7e1,g8a8a8dda67+a6fc98d2e7,g99cad8db69+f62e5b0af5,g9ddcbc5298+d4bad12328,ga1e77700b3+9c366c4306,ga8c6da7877+71e4819109,gb0e22166c9+25ba2f69a1,gb6a65358fc+469cd937ee,gbb8dafda3b+69d3c0e320,gc07e1c2157+a98bf949bb,gc120e1dc64+615ec43309,gc28159a63d+469cd937ee,gcf0d15dbbd+72a9f7b576,gdaeeff99f8+a38ce5ea23,ge6526c86ff+3a7c1ac5f1,ge79ae78c31+469cd937ee,gee10cc3b42+a6fc98d2e7,gf1cff7945b+21ad69e7e1,gfbcc870c63+9a11dc8c8f
LSST Data Management Base Package
Loading...
Searching...
No Matches
matchPessimisticB.py
Go to the documentation of this file.
1# This file is part of meas_astrom.
2#
3# Developed for the LSST Data Management System.
4# This product includes software developed by the LSST Project
5# (https://www.lsst.org).
6# See the COPYRIGHT file at the top-level directory of this distribution
7# for details of code ownership.
8#
9# This program is free software: you can redistribute it and/or modify
10# it under the terms of the GNU General Public License as published by
11# the Free Software Foundation, either version 3 of the License, or
12# (at your option) any later version.
13#
14# This program is distributed in the hope that it will be useful,
15# but WITHOUT ANY WARRANTY; without even the implied warranty of
16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17# GNU General Public License for more details.
18#
19# You should have received a copy of the GNU General Public License
20# along with this program. If not, see <https://www.gnu.org/licenses/>.
21
22__all__ = ["MatchPessimisticBTask", "MatchPessimisticBConfig",
23 "MatchTolerancePessimistic"]
24
25import numpy as np
26from scipy.spatial import cKDTree
27
28import lsst.pex.config as pexConfig
29import lsst.pipe.base as pipeBase
30import lsst.geom as geom
31import lsst.afw.table as afwTable
32from lsst.utils.timer import timeMethod
33
34from .matchOptimisticBTask import MatchTolerance
35
36from .pessimistic_pattern_matcher_b_3D import PessimisticPatternMatcherB
37
38
40 """Stores match tolerances for use in AstrometryTask and later
41 iterations of the matcher.
42
43 MatchPessimisticBTask relies on several state variables to be
44 preserved over different iterations in the
45 AstrometryTask.matchAndFitWcs loop of AstrometryTask.
46
47 Parameters
48 ----------
49 maxMatchDist : `lsst.geom.Angle`
50 Maximum distance to consider a match from the previous match/fit
51 iteration.
52 autoMaxMatchDist : `lsst.geom.Angle`
53 Automated estimation of the maxMatchDist from the sky statistics of the
54 source and reference catalogs.
55 maxShift : `lsst.geom.Angle`
56 Maximum shift found in the previous match/fit cycle.
57 lastMatchedPattern : `int`
58 Index of the last source pattern that was matched into the reference
59 data.
60 failedPatternList : `list` of `int`
61 Previous matches were found to be false positives.
62 PPMbObj : `lsst.meas.astrom.PessimisticPatternMatcherB`
63 Initialized Pessimistic pattern matcher object. Storing this prevents
64 the need for recalculation of the searchable distances in the PPMB.
65 """
66
67 def __init__(self, maxMatchDist=None, autoMaxMatchDist=None,
68 maxShift=None, lastMatchedPattern=None,
69 failedPatternList=None, PPMbObj=None):
70 self.maxMatchDistmaxMatchDist = maxMatchDist
71 self.autoMaxMatchDist = autoMaxMatchDist
72 self.maxShift = maxShift
73 self.lastMatchedPattern = lastMatchedPattern
74 self.PPMbObj = PPMbObj
75 if failedPatternList is None:
77 else:
78 self.failedPatternList = failedPatternList
79
80
81class MatchPessimisticBConfig(pexConfig.Config):
82 """Configuration for MatchPessimisticBTask
83 """
84 numBrightStars = pexConfig.RangeField(
85 doc="Maximum number of bright stars to use. Sets the max number of patterns "
86 "that can be tested.",
87 dtype=int,
88 default=150,
89 min=2,
90 )
91 minMatchedPairs = pexConfig.RangeField(
92 doc="Minimum number of matched pairs; see also minFracMatchedPairs.",
93 dtype=int,
94 default=30,
95 min=2,
96 )
97 minFracMatchedPairs = pexConfig.RangeField(
98 doc="Minimum number of matched pairs as a fraction of the smaller of "
99 "the number of reference stars or the number of good sources; "
100 "the actual minimum is the smaller of this value or "
101 "minMatchedPairs.",
102 dtype=float,
103 default=0.3,
104 min=0,
105 max=1,
106 )
107 matcherIterations = pexConfig.RangeField(
108 doc="Number of softening iterations in matcher.",
109 dtype=int,
110 default=5,
111 min=1,
112 )
113 maxOffsetPix = pexConfig.RangeField(
114 doc="Maximum allowed shift of WCS, due to matching (pixel). "
115 "When changing this value, the "
116 "LoadReferenceObjectsConfig.pixelMargin should also be updated.",
117 dtype=int,
118 default=250,
119 max=4000,
120 )
121 maxRotationDeg = pexConfig.RangeField(
122 doc="Rotation angle allowed between sources and position reference "
123 "objects (degrees).",
124 dtype=float,
125 default=1.0,
126 max=6.0,
127 )
128 numPointsForShape = pexConfig.Field(
129 doc="Number of points to define a shape for matching.",
130 dtype=int,
131 default=6,
132 )
133 numPointsForShapeAttempt = pexConfig.Field(
134 doc="Number of points to try for creating a shape. This value should "
135 "be greater than or equal to numPointsForShape. Besides "
136 "loosening the signal to noise cut in the 'matcher' SourceSelector, "
137 "increasing this number will solve CCDs where no match was found.",
138 dtype=int,
139 default=6,
140 )
141 minMatchDistPixels = pexConfig.RangeField(
142 doc="Distance in units of pixels to always consider a source-"
143 "reference pair a match. This prevents the astrometric fitter "
144 "from over-fitting and removing stars that should be matched and "
145 "allows for inclusion of new matches as the wcs improves.",
146 dtype=float,
147 default=1.0,
148 min=0.0,
149 max=6.0,
150 )
151 numPatternConsensus = pexConfig.Field(
152 doc="Number of implied shift/rotations from patterns that must agree "
153 "before it a given shift/rotation is accepted. This is only used "
154 "after the first softening iteration fails and if both the "
155 "number of reference and source objects is greater than "
156 "numBrightStars.",
157 dtype=int,
158 default=3,
159 )
160 numRefRequireConsensus = pexConfig.Field(
161 doc="If the available reference objects exceeds this number, "
162 "consensus/pessimistic mode will enforced regardless of the "
163 "number of available sources. Below this optimistic mode ("
164 "exit at first match rather than requiring numPatternConsensus to "
165 "be matched) can be used. If more sources are required to match, "
166 "decrease the signal to noise cut in the sourceSelector.",
167 dtype=int,
168 default=1000,
169 )
170 maxRefObjects = pexConfig.RangeField(
171 doc="Maximum number of reference objects to use for the matcher. The "
172 "absolute maximum allowed for is 2 ** 16 for memory reasons.",
173 dtype=int,
174 default=2**16,
175 min=0,
176 max=2**16 + 1,
177 )
178
179 def validate(self):
180 pexConfig.Config.validate(self)
182 raise ValueError("numPointsForShapeAttempt must be greater than "
183 "or equal to numPointsForShape.")
184 if self.numPointsForShape > self.numBrightStars:
185 raise ValueError("numBrightStars must be greater than "
186 "numPointsForShape.")
187
188
189# The following block adds links to this task from the Task Documentation page.
190# \addtogroup LSST_task_documentation
191# \{
192# \page measAstrom_MatchPessimisticBTask
193# \ref MatchPessimisticBTask "MatchPessimisticBTask"
194# Match sources to reference objects
195# \}
196
197
198class MatchPessimisticBTask(pipeBase.Task):
199 """Match sources to reference objects.
200 """
201
202 ConfigClass = MatchPessimisticBConfig
203 _DefaultName = "matchPessimisticB"
204
205 def __init__(self, **kwargs):
206 pipeBase.Task.__init__(self, **kwargs)
207
208 @timeMethod
209 def matchObjectsToSources(self, refCat, sourceCat, wcs, sourceFluxField, refFluxField,
210 match_tolerance=None):
211 """Match sources to position reference stars
212
213 refCat : `lsst.afw.table.SimpleCatalog`
214 catalog of reference objects that overlap the exposure; reads
215 fields for:
216
217 - coord
218 - the specified flux field
219
220 sourceCat : `lsst.afw.table.SourceCatalog`
221 Catalog of sources found on an exposure. This should already be
222 down-selected to "good"/"usable" sources in the calling Task.
223 wcs : `lsst.afw.geom.SkyWcs`
224 estimated WCS
225 sourceFluxField: `str`
226 field of sourceCat to use for flux
227 refFluxField : `str`
228 field of refCat to use for flux
229 match_tolerance : `lsst.meas.astrom.MatchTolerancePessimistic`
230 is a MatchTolerance class object or `None`. This this class is used
231 to communicate state between AstrometryTask and MatcherTask.
232 AstrometryTask will also set the MatchTolerance class variable
233 maxMatchDist based on the scatter AstrometryTask has found after
234 fitting for the wcs.
235
236 Returns
237 -------
238 result : `lsst.pipe.base.Struct`
239 Result struct with components:
240
241 - ``matches`` : source to reference matches found (`list` of
242 `lsst.afw.table.ReferenceMatch`)
243 - ``usableSourceCat`` : a catalog of sources potentially usable for
244 matching and WCS fitting (`lsst.afw.table.SourceCatalog`).
245 - ``match_tolerance`` : a MatchTolerance object containing the
246 resulting state variables from the match
247 (`lsst.meas.astrom.MatchTolerancePessimistic`).
248 """
249 import lsstDebug
250 debug = lsstDebug.Info(__name__)
251
252 # If we get an empty tolerance struct create the variables we need for
253 # this matcher.
254 if match_tolerance is None:
255 match_tolerance = MatchTolerancePessimistic()
256
257 # Make a name alias here for consistency with older code, and to make
258 # it clear that this is a good/usable (cleaned) source catalog.
259 goodSourceCat = sourceCat
260
261 numUsableSources = len(goodSourceCat)
262
263 if len(goodSourceCat) == 0:
264 raise pipeBase.TaskError("No sources are good")
265
266 minMatchedPairs = min(self.config.minMatchedPairs,
267 int(self.config.minFracMatchedPairs
268 * min([len(refCat), len(goodSourceCat)])))
269
270 if len(goodSourceCat) <= self.config.numPointsForShape:
271 msg = (f"Not enough catalog objects ({len(goodSourceCat)}) to make a "
272 f"shape for the matcher (need {self.config.numPointsForShape}).")
273 raise RuntimeError(msg)
274 if len(refCat) <= self.config.numPointsForShape:
275 msg = (f"Not enough refcat objects ({len(refCat)}) to make a "
276 f"shape for the matcher (need {self.config.numPointsForShape}).")
277 raise RuntimeError(msg)
278
279 if len(refCat) > self.config.maxRefObjects:
280 self.log.warning(
281 "WARNING: Reference catalog larger that maximum allowed. "
282 "Trimming to %i", self.config.maxRefObjects)
283 trimmedRefCat = self._filterRefCat(refCat, refFluxField)
284 else:
285 trimmedRefCat = refCat
286
287 doMatchReturn = self._doMatch(
288 refCat=trimmedRefCat,
289 sourceCat=goodSourceCat,
290 wcs=wcs,
291 refFluxField=refFluxField,
292 numUsableSources=numUsableSources,
293 minMatchedPairs=minMatchedPairs,
294 match_tolerance=match_tolerance,
295 sourceFluxField=sourceFluxField,
296 verbose=debug.verbose,
297 )
298 matches = doMatchReturn.matches
299 match_tolerance = doMatchReturn.match_tolerance
300
301 if len(matches) == 0:
302 raise RuntimeError("Unable to match sources")
303
304 self.log.info("Matched %d sources", len(matches))
305 if len(matches) < minMatchedPairs:
306 self.log.warning("Number of matches is smaller than request")
307
308 return pipeBase.Struct(
309 matches=matches,
310 usableSourceCat=goodSourceCat,
311 match_tolerance=match_tolerance,
312 )
313
314 def _filterRefCat(self, refCat, refFluxField):
315 """Sub-select a number of reference objects starting from the brightest
316 and maxing out at the number specified by maxRefObjects in the config.
317
318 No trimming is done if len(refCat) > config.maxRefObjects.
319
320 Parameters
321 ----------
322 refCat : `lsst.afw.table.SimpleCatalog`
323 Catalog of reference objects to trim.
324 refFluxField : `str`
325 field of refCat to use for flux
326 Returns
327 -------
328 outCat : `lsst.afw.table.SimpleCatalog`
329 Catalog trimmed to the number set in the task config from the
330 brightest flux down.
331 """
332 # Find the flux cut that gives us the desired number of objects.
333 if len(refCat) <= self.config.maxRefObjects:
334 return refCat
335 fluxArray = refCat.get(refFluxField)
336 sortedFluxArray = fluxArray[fluxArray.argsort()]
337 minFlux = sortedFluxArray[-(self.config.maxRefObjects + 1)]
338
339 selected = (refCat.get(refFluxField) > minFlux)
340
341 outCat = afwTable.SimpleCatalog(refCat.schema)
342 outCat.reserve(self.config.maxRefObjects)
343 outCat.extend(refCat[selected])
344
345 return outCat
346
347 @timeMethod
348 def _doMatch(self, refCat, sourceCat, wcs, refFluxField, numUsableSources,
349 minMatchedPairs, match_tolerance, sourceFluxField, verbose):
350 """Implementation of matching sources to position reference objects
351
352 Unlike matchObjectsToSources, this method does not check if the sources
353 are suitable.
354
355 Parameters
356 ----------
357 refCat : `lsst.afw.table.SimpleCatalog`
358 catalog of position reference objects that overlap an exposure
359 sourceCat : `lsst.afw.table.SourceCatalog`
360 catalog of sources found on the exposure
361 wcs : `lsst.afw.geom.SkyWcs`
362 estimated WCS of exposure
363 refFluxField : `str`
364 field of refCat to use for flux
365 numUsableSources : `int`
366 number of usable sources (sources with known centroid that are not
367 near the edge, but may be saturated)
368 minMatchedPairs : `int`
369 minimum number of matches
370 match_tolerance : `lsst.meas.astrom.MatchTolerancePessimistic`
371 a MatchTolerance object containing variables specifying matcher
372 tolerances and state from possible previous runs.
373 sourceFluxField : `str`
374 Name of the flux field in the source catalog.
375 verbose : `bool`
376 Set true to print diagnostic information to std::cout
377
378 Returns
379 -------
380 result :
381 Results struct with components:
382
383 - ``matches`` : a list the matches found
384 (`list` of `lsst.afw.table.ReferenceMatch`).
385 - ``match_tolerance`` : MatchTolerance containing updated values from
386 this fit iteration (`lsst.meas.astrom.MatchTolerancePessimistic`)
387 """
388
389 # Load the source and reference catalog as spherical points
390 # in numpy array. We do this rather than relying on internal
391 # lsst C objects for simplicity and because we require
392 # objects contiguous in memory. We need to do these slightly
393 # differently for the reference and source cats as they are
394 # different catalog objects with different fields.
395 src_array = np.empty((len(sourceCat), 4), dtype=np.float64)
396 for src_idx, srcObj in enumerate(sourceCat):
397 coord = wcs.pixelToSky(srcObj.getCentroid())
398 theta = np.pi / 2 - coord.getLatitude().asRadians()
399 phi = coord.getLongitude().asRadians()
400 flux = srcObj[sourceFluxField]
401 src_array[src_idx, :] = \
402 self._latlong_flux_to_xyz_mag(theta, phi, flux)
403
404 if match_tolerance.PPMbObj is None or \
405 match_tolerance.autoMaxMatchDist is None:
406 # The reference catalog is fixed per AstrometryTask so we only
407 # create the data needed if this is the first step in the match
408 # fit cycle.
409 ref_array = np.empty((len(refCat), 4), dtype=np.float64)
410 for ref_idx, refObj in enumerate(refCat):
411 theta = np.pi / 2 - refObj.getDec().asRadians()
412 phi = refObj.getRa().asRadians()
413 flux = refObj[refFluxField]
414 ref_array[ref_idx, :] = \
415 self._latlong_flux_to_xyz_mag(theta, phi, flux)
416 # Create our matcher object.
417 match_tolerance.PPMbObj = PessimisticPatternMatcherB(
418 ref_array[:, :3], self.log)
419 self.log.debug("Computing source statistics...")
420 maxMatchDistArcSecSrc = self._get_pair_pattern_statistics(
421 src_array)
422 self.log.debug("Computing reference statistics...")
423 maxMatchDistArcSecRef = self._get_pair_pattern_statistics(
424 ref_array)
425 maxMatchDistArcSec = np.max((
426 self.config.minMatchDistPixels
427 * wcs.getPixelScale().asArcseconds(),
428 np.min((maxMatchDistArcSecSrc,
429 maxMatchDistArcSecRef))))
430 match_tolerance.autoMaxMatchDist = geom.Angle(
431 maxMatchDistArcSec, geom.arcseconds)
432
433 # Set configurable defaults when we encounter None type or set
434 # state based on previous run of AstrometryTask._matchAndFitWcs.
435 if match_tolerance.maxShift is None:
436 maxShiftArcseconds = (self.config.maxOffsetPix
437 * wcs.getPixelScale().asArcseconds())
438 else:
439 # We don't want to clamp down too hard on the allowed shift so
440 # we test that the smallest we ever allow is the pixel scale.
441 maxShiftArcseconds = np.max(
442 (match_tolerance.maxShift.asArcseconds(),
443 self.config.minMatchDistPixels
444 * wcs.getPixelScale().asArcseconds()))
445
446 # If our tolerances are not set from a previous run, estimate a
447 # starting tolerance guess from the statistics of patterns we can
448 # create on both the source and reference catalog. We use the smaller
449 # of the two.
450 if match_tolerance.maxMatchDist is None:
451 match_tolerance.maxMatchDist = match_tolerance.autoMaxMatchDist
452 else:
453 maxMatchDistArcSec = np.max(
454 (self.config.minMatchDistPixels
455 * wcs.getPixelScale().asArcseconds(),
456 np.min((match_tolerance.maxMatchDist.asArcseconds(),
457 match_tolerance.autoMaxMatchDist.asArcseconds()))))
458
459 # Make sure the data we are considering is dense enough to require
460 # the consensus mode of the matcher. If not default to Optimistic
461 # pattern matcher behavior. We enforce pessimistic mode if the
462 # reference cat is sufficiently large, avoiding false positives.
463 numConsensus = self.config.numPatternConsensus
464 if len(refCat) < self.config.numRefRequireConsensus:
465 minObjectsForConsensus = \
466 self.config.numBrightStars + \
467 self.config.numPointsForShapeAttempt
468 if len(refCat) < minObjectsForConsensus or \
469 len(sourceCat) < minObjectsForConsensus:
470 numConsensus = 1
471
472 self.log.debug("Current tol maxDist: %.4f arcsec",
473 maxMatchDistArcSec)
474 self.log.debug("Current shift: %.4f arcsec",
475 maxShiftArcseconds)
476
477 match_found = False
478 # Start the iteration over our tolerances.
479 for soften_dist in range(self.config.matcherIterations):
480 if soften_dist == 0 and \
481 match_tolerance.lastMatchedPattern is not None:
482 # If we are on the first, most stringent tolerance,
483 # and have already found a match, the matcher should behave
484 # like an optimistic pattern matcher. Exiting at the first
485 # match.
486 run_n_consent = 1
487 else:
488 # If we fail or this is the first match attempt, set the
489 # pattern consensus to the specified config value.
490 run_n_consent = numConsensus
491 # We double the match dist tolerance each round and add 1 to the
492 # to the number of candidate spokes to check.
493 matcher_struct = match_tolerance.PPMbObj.match(
494 source_array=src_array,
495 n_check=self.config.numPointsForShapeAttempt,
496 n_match=self.config.numPointsForShape,
497 n_agree=run_n_consent,
498 max_n_patterns=self.config.numBrightStars,
499 max_shift=maxShiftArcseconds,
500 max_rotation=self.config.maxRotationDeg,
501 max_dist=maxMatchDistArcSec * 2. ** soften_dist,
502 min_matches=minMatchedPairs,
503 pattern_skip_array=np.array(
504 match_tolerance.failedPatternList)
505 )
506
507 if soften_dist == 0 and \
508 len(matcher_struct.match_ids) == 0 and \
509 match_tolerance.lastMatchedPattern is not None:
510 # If we found a pattern on a previous match-fit iteration and
511 # can't find an optimistic match on our first try with the
512 # tolerances as found in the previous match-fit,
513 # the match we found in the last iteration was likely bad. We
514 # append the bad match's index to the a list of
515 # patterns/matches to skip on subsequent iterations.
516 match_tolerance.failedPatternList.append(
517 match_tolerance.lastMatchedPattern)
518 match_tolerance.lastMatchedPattern = None
519 maxShiftArcseconds = \
520 self.config.maxOffsetPix * wcs.getPixelScale().asArcseconds()
521 elif len(matcher_struct.match_ids) > 0:
522 # Match found, save a bit a state regarding this pattern
523 # in the match tolerance class object and exit.
524 match_tolerance.maxShift = \
525 matcher_struct.shift * geom.arcseconds
526 match_tolerance.lastMatchedPattern = \
527 matcher_struct.pattern_idx
528 match_found = True
529 break
530
531 # If we didn't find a match, exit early.
532 if not match_found:
533 return pipeBase.Struct(
534 matches=[],
535 match_tolerance=match_tolerance,
536 )
537
538 # The matcher returns all the nearest neighbors that agree between
539 # the reference and source catalog. For the current astrometric solver
540 # we need to remove as many false positives as possible before sending
541 # the matches off to the solver. The low value of 100 and high value of
542 # 2 are the low number of sigma and high respectively. The exact values
543 # were found after testing on data of various reference/source
544 # densities and astrometric distortion quality, specifically the
545 # visits: HSC (3358), DECam (406285, 410827),
546 # CFHT (793169, 896070, 980526).
547 distances_arcsec = np.degrees(matcher_struct.distances_rad) * 3600
548 dist_cut_arcsec = np.max(
549 (np.degrees(matcher_struct.max_dist_rad) * 3600,
550 self.config.minMatchDistPixels * wcs.getPixelScale().asArcseconds()))
551
552 # A match has been found, return our list of matches and
553 # return.
554 matches = []
555 for match_id_pair, dist_arcsec in zip(matcher_struct.match_ids,
556 distances_arcsec):
557 if dist_arcsec < dist_cut_arcsec:
559 match.first = refCat[int(match_id_pair[1])]
560 match.second = sourceCat[int(match_id_pair[0])]
561 # We compute the true distance along and sphere. This isn't
562 # used in the WCS fitter however it is used in the unittest
563 # to confirm the matches computed.
564 match.distance = match.first.getCoord().separation(
565 match.second.getCoord()).asArcseconds()
566 matches.append(match)
567
568 return pipeBase.Struct(
569 matches=matches,
570 match_tolerance=match_tolerance,
571 )
572
573 def _latlong_flux_to_xyz_mag(self, theta, phi, flux):
574 """Convert angles theta and phi and a flux into unit sphere
575 x, y, z, and a relative magnitude.
576
577 Takes in a afw catalog object and converts the catalog object RA, DECs
578 to points on the unit sphere. Also converts the flux into a simple,
579 non-zero-pointed magnitude for relative sorting.
580
581 Parameters
582 ----------
583 theta : `float`
584 Angle from the north pole (z axis) of the sphere
585 phi : `float`
586 Rotation around the sphere
587
588 Return
589 ------
590 output_array : `numpy.ndarray`, (N, 4)
591 Spherical unit vector x, y, z with flux.
592 """
593 output_array = np.empty(4, dtype=np.float64)
594 output_array[0] = np.sin(theta)*np.cos(phi)
595 output_array[1] = np.sin(theta)*np.sin(phi)
596 output_array[2] = np.cos(theta)
597 if flux > 0:
598 output_array[3] = -2.5 * np.log10(flux)
599 else:
600 # Set flux to a very faint mag if its for some reason it
601 # does not exist
602 output_array[3] = 99.
603
604 return output_array
605
606 def _get_pair_pattern_statistics(self, cat_array):
607 """ Compute the tolerances for the matcher automatically by comparing
608 pinwheel patterns as we would in the matcher.
609
610 We test how similar the patterns we can create from a given set of
611 objects by computing the spoke lengths for each pattern and sorting
612 them from smallest to largest. The match tolerance is the average
613 distance per spoke between the closest two patterns in the sorted
614 spoke length space.
615
616 Parameters
617 ----------
618 cat_array : `numpy.ndarray`, (N, 3)
619 array of 3 vectors representing the x, y, z position of catalog
620 objects on the unit sphere.
621
622 Returns
623 -------
624 dist_tol : `float`
625 Suggested max match tolerance distance calculated from comparisons
626 between pinwheel patterns used in optimistic/pessimistic pattern
627 matcher.
628 """
629
630 self.log.debug("Starting automated tolerance calculation...")
631
632 # Create an empty array of all the patterns we possibly make
633 # sorting from brightest to faintest.
634 pattern_array = np.empty(
635 (cat_array.shape[0] - self.config.numPointsForShape,
636 self.config.numPointsForShape - 1))
637 flux_args_array = np.argsort(cat_array[:, -1])
638
639 # Sort our input array.
640 tmp_sort_array = cat_array[flux_args_array]
641
642 # Start making patterns.
643 for start_idx in range(cat_array.shape[0]
644 - self.config.numPointsForShape):
645 pattern_points = tmp_sort_array[start_idx:start_idx
646 + self.config.numPointsForShape, :-1]
647 pattern_delta = pattern_points[1:, :] - pattern_points[0, :]
648 pattern_array[start_idx, :] = np.sqrt(
649 pattern_delta[:, 0] ** 2
650 + pattern_delta[:, 1] ** 2
651 + pattern_delta[:, 2] ** 2)
652
653 # When we store the length of each spoke in our pattern we
654 # sort from shortest to longest so we have a defined space
655 # to compare them in.
656 pattern_array[start_idx, :] = pattern_array[
657 start_idx, np.argsort(pattern_array[start_idx, :])]
658
659 # Create a searchable tree object of the patterns and find
660 # for any given pattern the closest pattern in the sorted
661 # spoke length space.
662 dist_tree = cKDTree(
663 pattern_array[:, :(self.config.numPointsForShape - 1)])
664 dist_nearest_array, ids = dist_tree.query(
665 pattern_array[:, :(self.config.numPointsForShape - 1)], k=2)
666 dist_nearest_array = dist_nearest_array[:, 1]
667 dist_nearest_array.sort()
668
669 # We use the two closest patterns to set our tolerance.
670 dist_idx = 0
671 dist_tol = (np.degrees(dist_nearest_array[dist_idx]) * 3600.
672 / (self.config.numPointsForShape - 1.))
673
674 self.log.debug("Automated tolerance")
675 self.log.debug("\tdistance/match tol: %.4f [arcsec]", dist_tol)
676
677 return dist_tol
int min
Tag types used to declare specialized field types.
Definition misc.h:31
Custom catalog class for record/table subclasses that are guaranteed to have an ID,...
A class representing an angle.
Definition Angle.h:128
_doMatch(self, refCat, sourceCat, wcs, refFluxField, numUsableSources, minMatchedPairs, match_tolerance, sourceFluxField, verbose)
matchObjectsToSources(self, refCat, sourceCat, wcs, sourceFluxField, refFluxField, match_tolerance=None)
__init__(self, maxMatchDist=None, autoMaxMatchDist=None, maxShift=None, lastMatchedPattern=None, failedPatternList=None, PPMbObj=None)