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