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