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