lsst.meas.astrom  15.0-8-gb6a5594
pessimistic_pattern_matcher_b_3D.py
Go to the documentation of this file.
1 
2 
3 import numpy as np
4 from scipy.optimize import least_squares
5 from scipy.spatial import cKDTree
6 
7 import lsst.pipe.base as pipeBase
8 
9 
10 def _rotation_matrix_chi_sq(flattened_rot_matrix,
11  pattern_a,
12  pattern_b,
13  max_dist_rad):
14  """Compute the squared differences for least squares fitting.
15 
16  Given a flattened rotation matrix, one N point pattern and another N point
17  pattern to transform into to, compute the squared differences between the
18  points in the two patterns after the rotation.
19 
20  Parameters
21  ----------
22  flattened_rot_matrix : float array
23  A flattened array representing a 3x3 rotation matrix. The array is
24  flattened to comply with the API of scipy.optimize.least_squares.
25  Flattened elements are [[0, 0], [0, 1], [0, 2], [1, 0]...]
26  pattern_a : float array of 3 vectors
27  A array containing N, 3 vectors representing the objects we would like
28  to transform into the frame of pattern_b.
29  pattern_b : float array of 3 vectors
30  A array containing N, 3 vectors representing the reference frame we
31  would like to transform pattern_a into.
32  max_dist_rad : float
33  The maximum distance allowed from the pattern matching. This value is
34  used as the standard error for the resultant chi values.
35 
36  Returns
37  -------
38  float array
39  Array of differences between the vectors representing of the source
40  pattern rotated into the reference frame and the converse. This is
41  used to minimize in a least squares fitter.
42  """
43  # Unflatten the rotation matrix
44  rot_matrix = flattened_rot_matrix.reshape((3, 3))
45  # Compare the rotated source pattern to the references.
46  rot_pattern_a = np.dot(rot_matrix, pattern_a.transpose()).transpose()
47  diff_pattern_a_to_b = rot_pattern_a - pattern_b
48  # Compute the matrix inverse and compare the opposite transform to better
49  # constrain the matrix.
50  rot_pattern_b = np.dot(np.linalg.inv(rot_matrix),
51  pattern_b.transpose()).transpose()
52  diff_pattern_b_to_a = rot_pattern_b - pattern_a
53  # Return the flattened differences and length tolerances for use in a least
54  # squares fitter.
55  return np.concatenate(
56  (diff_pattern_a_to_b.flatten() / max_dist_rad,
57  diff_pattern_b_to_a.flatten() / max_dist_rad))
58 
59 
61  """ Class implementing a pessimistic version of Optimistic Pattern Matcher
62  B (OPMb) from Tabur 2007. The class loads and stores the reference object
63  in a convenient data structure for matching any set of source objects that
64  are assumed to contain each other. The pessimistic nature of the algorithm
65  comes from requiring that it discovers at least two patterns that agree on
66  the correct shift and rotation for matching before exiting. The original
67  behavior of OPMb can be recovered simply. Patterns matched between the
68  input datasets are n-spoked pinwheels created from n+1 points. Refer to
69  DMTN #031 for more details. http://github.com/lsst-dm/dmtn-031
70  --------------------------------------------------------------------------
71  Attributes:
72  reference_array : float array
73  spherical points x, y, z of to use as reference objects for
74  pattern matching.
75  log : an lsst.log instance
76  pair_id_array : int array
77  Internal lookup table. Given an id in the reference array, return
78  an array of the id pair that contains this object's id sorted on
79  the distance to the pairs.
80  pair_delta_array : float array
81  Internal lookup table. Given an id in the reference array, return
82  an array of the 3 vector deltas of all other pairs sorted on their
83  distance.
84  pair_dist_array : float array
85  Internal lookup table. Given an id in the reference return an array
86  of pair distances of all other pairs sorted on distance.
87  dist_array : float array
88  Array of all pairs of objects in the reference array sorted on
89  distance.
90  id_array : int array
91  Array of id pairs that lookup into the reference array sorted
92  on pair distance.
93  delta_array : float array
94  Array of 3 vector deltas for each pair in the reference array
95  sorted on pair distance.
96  """
97 
98  def __init__(self, reference_array, log):
99  """
100  Arguments
101  ---------
102  reference_array : float array
103  Array of spherical points x, y, z to use as reference objects.
104  log : lsst.log
105  logger object for reporting warnings and failures.
106  """
107  self._reference_array = reference_array
108  self._n_reference = len(self._reference_array)
109  self.log = log
110 
112 
113  def _build_distances_and_angles(self):
114  """ Create the data structures we will use to search for our pattern
115  match in.
116 
117  Throughout this function and the rest of the
118 
119  class we use id to reference the position in the input reference
120  catalog and index to 'index' into the arrays sorted on distance.
121 
122  Parameters
123  ----------
124  self
125 
126  Returns
127  -------
128  None
129  """
130 
131  # Initialize the arrays we will need for quick look up of pairs once
132  # have a candidate spoke center.
133  self._pair_id_array = np.empty(
134  (self._n_reference, self._n_reference - 1),
135  dtype=np.uint32)
136  self._pair_delta_array = np.empty(
137  (self._n_reference, self._n_reference - 1, 3),
138  dtype=np.float64)
139  self._pair_dist_array = np.empty(
140  (self._n_reference, self._n_reference - 1),
141  dtype=np.float64)
142 
143  # Create empty lists to temporarily store our pair information per
144  # reference object. These will be concatenated into our final arrays.
145  sub_id_array_list = []
146  sub_delta_array_list = []
147  sub_dist_array_list = []
148 
149  # Loop over reference objects and store pair distances, ids, and
150  # 3 vector deltas.
151  for ref_id, ref_obj in enumerate(self._reference_array):
152 
153  # Reserve and fill the ids of each reference object pair.
154  sub_id_array = np.zeros((self._n_reference - 1 - ref_id, 2),
155  dtype=np.uint32)
156  sub_id_array[:, 0] = ref_id
157  sub_id_array[:, 1] = np.arange(ref_id + 1, self._n_reference,
158  dtype=np.uint32)
159 
160  # Compute the vector deltas for each pair of reference objects
161  # and compute and store the distances.
162  sub_delta_array = (self._reference_array[ref_id + 1:, :] -
163  ref_obj)
164  sub_dist_array = np.sqrt(sub_delta_array[:, 0] ** 2 +
165  sub_delta_array[:, 1] ** 2 +
166  sub_delta_array[:, 2] ** 2)
167 
168  # Append to our arrays to the output lists for later
169  # concatenation.
170  sub_id_array_list.append(sub_id_array)
171  sub_delta_array_list.append(sub_delta_array)
172  sub_dist_array_list.append(sub_dist_array)
173 
174  # Fill the pair look up arrays row wise and then column wise.
175  self._pair_id_array[ref_id, ref_id:] = sub_id_array[:, 1]
176  self._pair_delta_array[ref_id, ref_id:, :] = sub_delta_array
177  self._pair_dist_array[ref_id, ref_id:] = sub_dist_array
178 
179  # Don't fill the array column wise if we are on the last as
180  # these pairs have already been filled by previous
181  # iterations.
182  if ref_id < self._n_reference - 1:
183  self._pair_id_array[ref_id + 1:, ref_id] = ref_id
184  self._pair_delta_array[ref_id + 1:, ref_id, :] = \
185  sub_delta_array
186  self._pair_dist_array[ref_id + 1:, ref_id] = sub_dist_array
187 
188  # Sort each row on distance for fast look up of pairs given
189  # the id of one of the objects in the pair.
190  sorted_pair_dist_args = self._pair_dist_array[ref_id, :].argsort()
191  self._pair_dist_array[ref_id, :] = self._pair_dist_array[
192  ref_id, sorted_pair_dist_args]
193  self._pair_id_array[ref_id, :] = self._pair_id_array[
194  ref_id, sorted_pair_dist_args]
195  self._pair_delta_array[ref_id, :, :] = self._pair_delta_array[
196  ref_id, sorted_pair_dist_args, :]
197 
198  # Concatenate our arrays together.
199  unsorted_id_array = np.concatenate(sub_id_array_list)
200  unsorted_delta_array = np.concatenate(sub_delta_array_list)
201  unsorted_dist_array = np.concatenate(sub_dist_array_list)
202 
203  # Sort each array on the pair distances for the initial
204  # optimistic pattern matcher lookup.
205  sorted_dist_args = unsorted_dist_array.argsort()
206  self._dist_array = unsorted_dist_array[sorted_dist_args]
207  self._id_array = unsorted_id_array[sorted_dist_args]
208  self._delta_array = unsorted_delta_array[sorted_dist_args]
209 
210  return None
211 
212  def match(self, source_array, n_check, n_match, n_agree,
213  max_n_patterns, max_shift, max_rotation, max_dist,
214  min_matches, pattern_skip_array=None):
215  """Match a given source catalog into the loaded reference catalog.
216 
217  Given array of points on the unit sphere and tolerances, we
218  attempt to match a pinwheel like pattern between these input sources
219  and the reference objects this class was created with. This pattern
220  informs of the shift and rotation needed to align the input source
221  objects into the frame of the references.
222 
223  Parameters
224  ----------
225  source_array: float array
226  An array of spherical x,y,z coordinates and a magnitude in units
227  of objects having a lower value for sorting. The array should be
228  of shape (N, 4).
229  n_check : int value
230  Number of sources to create a pattern from. Not all objects may be
231  checked if n_match criteria is before looping through all n_check
232  objects.
233  n_match : int value
234  Number of objects to use in constructing a pattern to match.
235  n_agree: int value
236  Number of found patterns that must agree on their shift and
237  rotation before exiting. Set this value to 1 to recover the
238  expected behavior of Optimistic Pattern Matcher B.
239  max_n_patters : int value
240  Number of patterns to create from the input source objects to
241  attempt to match into the reference objects.
242  max_shift: float value
243  Maximum allowed shift to match patterns in arcseconds.
244  max_rotation: float value
245  Maximum allowed rotation between patterns in degrees.
246  max_dist: float value
247  Maximum distance in arcseconds allowed between candidate spokes in
248  the source and reference objects. Also sets that maximum distance
249  in the intermediate verify, pattern shift/rotation agreement, and
250  final verify steps.
251  pattern_skip_array: int array
252  Patterns we would like to skip. This could be due to the pattern
253  being matched on a previous iteration that we now consider invalid.
254  This assumes the ordering of the source objects is the same
255  between different runs of the matcher which, assuming no object
256  has been inserted or the magnitudes have changed, it should be.
257 
258  Returns
259  -------
260  output_struct : pipe.base.struct
261  A lsst.pipe.base struct containing the following outputs.
262 
263  matches : int array
264  (N, 2) array of matched ids for pairs. Empty list if no
265  match found.
266  distances_rad : float array
267  Radian distances between the matched objects. Empty list
268  if no match found.
269  pattern_idx : int
270  Index of matched pattern. None if no match found.
271  shift : float
272  Magnitude for the shift between the source and
273  reference objects in arcseconds. None if no match found.
274  """
275 
276  # Given our input source_array we sort on magnitude.
277  sorted_source_array = source_array[source_array[:, -1].argsort(), :3]
278  n_source = len(sorted_source_array)
279 
280  # Initialize output struct.
281  output_match_struct = pipeBase.Struct(
282  match_ids=[],
283  distances_rad=[],
284  pattern_idx=None,
285  shift=None,)
286 
287  if n_source <= 0:
288  self.log.warn("Source object array is empty. Unable to match. "
289  "Exiting matcher.")
290  return None
291 
292  # To test if the shifts and rotations we find agree with each other,
293  # we first create two test points situated at the top and bottom of
294  # where the z axis on the sphere bisects the source catalog.
295  test_vectors = self._compute_test_vectors(source_array[:, :3])
296 
297  # We now create an empty list of our resultant rotated vectors to
298  # compare the different rotations we find.
299  rot_vect_list = []
300 
301  # Convert the tolerances to values we will use in the code.
302  max_cos_shift = np.cos(np.radians(max_shift / 3600.))
303  max_cos_rot_sq = np.cos(np.radians(max_rotation)) ** 2
304  max_dist_rad = np.radians(max_dist / 3600.)
305 
306  # Loop through the sources from brightest to faintest, grabbing a
307  # chunk of n_check each time.
308  for pattern_idx in range(np.min((max_n_patterns,
309  n_source - n_match))):
310 
311  # If this pattern is one that we matched on the past but we
312  # now want to skip, we do so here.
313  if pattern_skip_array is not None and \
314  np.any(pattern_skip_array == pattern_idx):
315  self.log.debug(
316  "Skipping previously matched bad pattern %i..." %
317  pattern_idx)
318  continue
319  # Grab the sources to attempt to create this pattern.
320  pattern = sorted_source_array[
321  pattern_idx: np.min((pattern_idx + n_check, n_source)), :3]
322 
323  # Construct a pattern given the number of points defining the
324  # pattern complexity. This is the start of the primary tests to
325  # match our source pattern into the reference objects.
326  construct_return_struct = \
328  pattern, n_match, max_cos_shift, max_cos_rot_sq,
329  max_dist_rad)
330 
331  # Our struct is None if we could not match the pattern.
332  if construct_return_struct.ref_candidates is None or \
333  construct_return_struct.shift_rot_matrix is None or \
334  construct_return_struct.cos_shift is None or \
335  construct_return_struct.sin_rot is None:
336  continue
337 
338  # Grab the output data from the Struct object.
339  ref_candidates = construct_return_struct.ref_candidates
340  shift_rot_matrix = construct_return_struct.shift_rot_matrix
341  cos_shift = construct_return_struct.cos_shift
342  sin_rot = construct_return_struct.sin_rot
343 
344  # If we didn't match enough candidates we continue to the next
345  # pattern.
346  if len(ref_candidates) < n_match:
347  continue
348 
349  # Now that we know our pattern and shift/rotation are valid we
350  # store the the rotated versions of our test points for later
351  # use.
352  tmp_rot_vect_list = []
353  for test_vect in test_vectors:
354  tmp_rot_vect_list.append(np.dot(shift_rot_matrix, test_vect))
355  # Since our test point are in the source frame, we can test if
356  # their lengths are mostly preserved under the transform.
357  if not self._test_pattern_lengths(np.array(tmp_rot_vect_list),
358  max_dist_rad):
359  continue
360 
361  tmp_rot_vect_list.append(pattern_idx)
362  rot_vect_list.append(tmp_rot_vect_list)
363 
364  # Test if we have enough rotations, which agree, or if we
365  # are in optimistic mode.
366  if self._test_rotation_agreement(rot_vect_list, max_dist_rad) < \
367  n_agree - 1:
368  continue
369 
370  # Perform final verify.
371  match_sources_struct = self._match_sources(source_array[:, :3],
372  shift_rot_matrix)
373 
374  n_matched = len(match_sources_struct.match_ids[
375  match_sources_struct.distances_rad < max_dist_rad])
376 
377  # Check that we have enough matches.
378  if n_matched >= min_matches:
379  # Convert the observed shift to arcseconds
380  shift = np.degrees(np.arccos(cos_shift)) * 3600.
381  # Print information to the logger.
382  self.log.debug("Succeeded after %i patterns." % pattern_idx)
383  self.log.debug("\tShift %.4f arcsec" % shift)
384  self.log.debug("\tRotation: %.4f deg" %
385  np.degrees(np.arcsin(sin_rot)))
386 
387  # Fill the struct and return.
388  output_match_struct.match_ids = \
389  match_sources_struct.match_ids
390  output_match_struct.distances_rad = \
391  match_sources_struct.distances_rad
392  output_match_struct.pattern_idx = pattern_idx
393  output_match_struct.shift = shift
394  return output_match_struct
395 
396  self.log.warn("Failed after %i patterns." % (pattern_idx + 1))
397  return output_match_struct
398 
399  def _compute_test_vectors(self, source_array):
400  """Compute spherical 3 vectors at the edges of the x, y, z extent
401  of the input source catalog.
402 
403  Parameters
404  ----------
405  source_array : float array (N, 3)
406  array of 3 vectors representing positions on the unit
407  sphere.
408 
409  Returns
410  -------
411  float array of 3 vectors
412  Array of vectors representing the maximum extents in x, y, z
413  of the input source array. These are used with the rotations
414  the code finds to test for agreement from different patterns
415  when the code is running in pessimistic mode.
416  """
417 
418  # Get the center of source_array.
419  if np.any(np.logical_not(np.isfinite(source_array))):
420  self.log.warn("Input source objects contain non-finite values. "
421  "This could end badly.")
422  center_vect = np.nanmean(source_array, axis=0)
423 
424  # So that our rotation test works over the full sky we compute
425  # the max extent in each Cartesian direction x,y,z.
426  xbtm_vect = np.array([np.min(source_array[:, 0]), center_vect[1],
427  center_vect[2]], dtype=np.float64)
428  xtop_vect = np.array([np.max(source_array[:, 0]), center_vect[1],
429  center_vect[2]], dtype=np.float64)
430  xbtm_vect /= np.sqrt(np.dot(xbtm_vect, xbtm_vect))
431  xtop_vect /= np.sqrt(np.dot(xtop_vect, xtop_vect))
432 
433  ybtm_vect = np.array([center_vect[0], np.min(source_array[:, 1]),
434  center_vect[2]], dtype=np.float64)
435  ytop_vect = np.array([center_vect[0], np.max(source_array[:, 1]),
436  center_vect[2]], dtype=np.float64)
437  ybtm_vect /= np.sqrt(np.dot(ybtm_vect, ybtm_vect))
438  ytop_vect /= np.sqrt(np.dot(ytop_vect, ytop_vect))
439 
440  zbtm_vect = np.array([center_vect[0], center_vect[1],
441  np.min(source_array[:, 2])], dtype=np.float64)
442  ztop_vect = np.array([center_vect[0], center_vect[1],
443  np.max(source_array[:, 2])], dtype=np.float64)
444  zbtm_vect /= np.sqrt(np.dot(zbtm_vect, zbtm_vect))
445  ztop_vect /= np.sqrt(np.dot(ztop_vect, ztop_vect))
446 
447  # Return our list of vectors for later rotation testing.
448  return np.array([xbtm_vect, xtop_vect, ybtm_vect, ytop_vect,
449  zbtm_vect, ztop_vect])
450 
451  def _construct_pattern_and_shift_rot_matrix(self, src_pattern_array,
452  n_match, max_cos_theta_shift,
453  max_cos_rot_sq, max_dist_rad):
454  """Test an input source pattern against the reference catalog.
455 
456  Returns the candidate matched patterns and their
457  implied rotation matrices or None.
458 
459  Parameters
460  ----------
461  src_pattern_array : float array
462  Sub selection of source 3 vectors to create a pattern from
463  n_match : int
464  Number of points to attempt to create a pattern from. Must be
465  >= len(src_pattern_array)
466  max_cos_theta_shift : float
467  Maximum shift allowed between two patterns' centers.
468  max_cos_rot_sq : float
469  Maximum rotation between two patterns that have been shifted
470  to have their centers on top of each other.
471  max_dist_rad : float
472  Maximum delta distance allowed between the source and reference
473  pair distances to consider the reference pair a candidate for
474  the source pair. Also sets the tolerance between the opening
475  angles of the spokes when compared to the reference.
476 
477  Return
478  -------
479  lsst.pipe.base.Struct
480  Return a Struct containing the following data:
481 
482  ref_candidates : list of ints
483  ids of the matched pattern in the internal reference_array
484  object.
485  src_candidates : list of ints
486  Pattern ids of the sources matched.
487  shift_rot_matrix_src_to_ref : float array
488  3x3 matrix specifying the full shift and rotation between the
489  reference and source objects. Rotates
490  source into reference frame. None if match is not found.
491  shift_rot_matrix_ref_to_src : float array
492  3x3 matrix specifying the full shift
493  and rotation of the reference and source objects. Rotates
494  reference into source frame. None if match is not found.
495  cos_shift : float
496  Magnitude of the shift found between the two patten
497  centers. None if match is not found.
498  sin_rot : float value of the rotation to align the already shifted
499  source pattern to the reference pattern. None if no match
500  found.
501  """
502 
503  # Create our place holder variables for the matched sources and
504  # references. The source list starts with the 0th and first indexed
505  # objects as we are guaranteed to use those and these define both
506  # the shift and rotation of the final pattern.
507  output_matched_pattern = pipeBase.Struct(
508  ref_candidates=[],
509  src_candidates=[],
510  shift_rot_matrix=None,
511  cos_shift=None,
512  sin_rot=None)
513 
514  # Create the delta vectors and distances we will need to assemble the
515  # spokes of the pattern.
516  src_delta_array = np.empty((len(src_pattern_array) - 1, 3))
517  src_delta_array[:, 0] = (src_pattern_array[1:, 0] -
518  src_pattern_array[0, 0])
519  src_delta_array[:, 1] = (src_pattern_array[1:, 1] -
520  src_pattern_array[0, 1])
521  src_delta_array[:, 2] = (src_pattern_array[1:, 2] -
522  src_pattern_array[0, 2])
523  src_dist_array = np.sqrt(src_delta_array[:, 0]**2 +
524  src_delta_array[:, 1]**2 +
525  src_delta_array[:, 2]**2)
526 
527  # Our first test. We search the reference dataset for pairs
528  # that have the same length as our first source pairs to with
529  # plus/minus the max_dist tolerance.
530  ref_dist_index_array = self._find_candidate_reference_pairs(
531  src_dist_array[0], self._dist_array, max_dist_rad)
532 
533  # Start our loop over the candidate reference objects.
534  for ref_dist_idx in ref_dist_index_array:
535  # We have two candidates for which reference object corresponds
536  # with the source at the center of our pattern. As such we loop
537  # over and test both possibilities.
538  tmp_ref_pair_list = self._id_array[ref_dist_idx]
539  for pair_idx, ref_id in enumerate(tmp_ref_pair_list):
540  src_candidates = [0, 1]
541  ref_candidates = []
542  shift_rot_matrix = None
543  cos_shift = None
544  sin_rot = None
545  # Test the angle between our candidate ref center and the
546  # source center of our pattern. This angular distance also
547  # defines the shift we will later use.
548  ref_center = self._reference_array[ref_id]
549  cos_shift = np.dot(src_pattern_array[0], ref_center)
550  if cos_shift < max_cos_theta_shift:
551  continue
552 
553  # We can now append this one as a candidate.
554  ref_candidates.append(ref_id)
555  ref_delta = self._delta_array[ref_dist_idx]
556  # If the candidate reference center we found is second in
557  # this pair we need to reverse the direction of the
558  # corresponding pair's delta vector.
559  if pair_idx == 0:
560  ref_candidates.append(
561  tmp_ref_pair_list[1])
562  else:
563  ref_candidates.append(
564  tmp_ref_pair_list[0])
565  ref_delta *= -1
566 
567  # For dense fields it will be faster to compute the absolute
568  # rotation this pair suggests first rather than saving it
569  # after all the spokes are found. We then compute the cos^2
570  # of the rotation and first part of the rotation matrix from
571  # source to reference frame.
572  test_rot_struct = self._test_rotation(
573  src_pattern_array[0], ref_center, src_delta_array[0],
574  ref_delta, cos_shift, max_cos_rot_sq)
575  if test_rot_struct.cos_rot_sq is None or \
576  test_rot_struct.shift_matrix is None:
577  continue
578 
579  # Get the data from the return struct.
580  cos_rot_sq = test_rot_struct.cos_rot_sq
581  shift_matrix = test_rot_struct.shift_matrix
582 
583  # Now that we have a candidate first spoke and reference
584  # pattern center, we mask our future search to only those
585  # pairs that contain our candidate reference center.
586  ref_dist = self._dist_array[ref_dist_idx]
587  tmp_ref_delta_array = self._pair_delta_array[ref_id]
588  tmp_ref_dist_arary = self._pair_dist_array[ref_id]
589  tmp_ref_id_array = self._pair_id_array[ref_id]
590 
591  # Now we feed this sub data to match the spokes of
592  # our pattern.
593  pattern_spoke_struct = self._create_pattern_spokes(
594  src_pattern_array[0], src_delta_array, src_dist_array,
595  self._reference_array[ref_id], ref_id, ref_delta,
596  ref_dist, tmp_ref_delta_array, tmp_ref_dist_arary,
597  tmp_ref_id_array, max_dist_rad,
598  n_match)
599 
600  # If we don't find enough candidates we can continue to the
601  # next reference center pair.
602  if len(pattern_spoke_struct.ref_spoke_list) < n_match - 2 or \
603  len(pattern_spoke_struct.src_spoke_list) < n_match - 2:
604  continue
605 
606  # If we have the right number of matched ids we store these.
607  ref_candidates.extend(pattern_spoke_struct.ref_spoke_list)
608  src_candidates.extend(pattern_spoke_struct.src_spoke_list)
609 
610  # We can now create our full rotation matrix for both the
611  # shift and rotation. Reminder shift, aligns the pattern
612  # centers, rotation rotates the spokes on top of each other.
613  shift_rot_struct = self._create_shift_rot_matrix(
614  cos_rot_sq, shift_matrix, src_delta_array[0],
615  self._reference_array[ref_id], ref_delta)
616  # If we fail to create the rotation matrix, continue to the
617  # next objects.
618  if shift_rot_struct.sin_rot is None or \
619  shift_rot_struct.shift_rot_matrix is None:
620  continue
621 
622  # Get the data from the return struct.
623  sin_rot = shift_rot_struct.sin_rot
624  shift_rot_matrix = shift_rot_struct.shift_rot_matrix
625 
626  # Now that we have enough candidates we test to see if it
627  # passes intermediate verify. This shifts and rotates the
628  # source pattern into the reference frame and then tests that
629  # each source/reference object pair is within max_dist. It also
630  # tests the opposite rotation that is reference to source
631  # frame.
632  fit_shift_rot_matrix = self._intermediate_verify(
633  src_pattern_array[src_candidates],
634  self._reference_array[ref_candidates],
635  shift_rot_matrix, max_dist_rad)
636 
637  if fit_shift_rot_matrix is not None:
638  # Fill the struct and return.
639  output_matched_pattern.ref_candidates = ref_candidates
640  output_matched_pattern.src_candidates = src_candidates
641  output_matched_pattern.shift_rot_matrix = \
642  fit_shift_rot_matrix
643  output_matched_pattern.cos_shift = cos_shift
644  output_matched_pattern.sin_rot = sin_rot
645  return output_matched_pattern
646 
647  return output_matched_pattern
648 
649  def _find_candidate_reference_pairs(self, src_dist, ref_dist_array,
650  max_dist_rad):
651  """Wrap numpy.searchsorted to find the range of reference spokes
652  within a spoke distance tolerance of our source spoke.
653 
654  Returns an array sorted from the smallest absolute delta distance
655  between source and reference spoke length. This sorting increases the
656  speed for the pattern search greatly.
657 
658  Parameters
659  ----------
660  src_dist : float radians
661  float value of the distance we would like to search for in
662  the reference array in radians.
663  ref_dist_array : float array
664  sorted array of distances in radians.
665  max_dist_rad : float
666  maximum plus/minus search to find in the reference array in
667  radians.
668 
669  Return
670  ------
671  int array
672  indices lookup into the input ref_dist_array sorted by the
673  difference in value to the src_dist from absolute value
674  smallest to largest.
675  """
676  # Find the index of the minimum and maximum values that satisfy
677  # the tolerance.
678  start_idx = np.searchsorted(ref_dist_array, src_dist - max_dist_rad)
679  end_idx = np.searchsorted(ref_dist_array, src_dist + max_dist_rad,
680  side='right')
681 
682  # If these are equal there are no candidates and we exit.
683  if start_idx == end_idx:
684  return []
685 
686  # Make sure the endpoints of the input array are respected.
687  if start_idx < 0:
688  start_idx = 0
689  if end_idx > ref_dist_array.shape[0]:
690  end_idx = ref_dist_array.shape[0]
691 
692  # Now we sort the indices from smallest absolute delta dist difference
693  # to the largest and return the vector. This step greatly increases
694  # the speed of the algorithm.
695  tmp_diff_array = np.fabs(ref_dist_array[start_idx:end_idx] - src_dist)
696  return tmp_diff_array.argsort() + start_idx
697 
698  def _test_rotation(self, src_center, ref_center, src_delta, ref_delta,
699  cos_shift, max_cos_rot_sq):
700  """ Test if the rotation implied between the source
701  pattern and reference pattern is within tolerance. To test this
702  we need to create the first part of our spherical rotation matrix
703  which we also return for use later.
704 
705  Parameters
706  ----------
707  src_center : float array3
708  pattern.
709  ref_center : float array
710  3 vector defining the center of the candidate reference pinwheel
711  pattern.
712  src_delta : float array
713  3 vector delta between the source pattern center and the end of
714  the pinwheel spoke.
715  ref_delta : float array
716  3 vector delta of the candidate matched reference pair
717  cos_shift : float
718  Cosine of the angle between the source and reference candidate
719  centers.
720  max_cos_rot_sq : float
721  candidate reference pair after shifting the centers on top of each
722  other. The function will return None if the rotation implied is
723  greater than max_cos_rot_sq.
724 
725  Returns
726  -------
727  lsst.pipe.base.Struct
728  Return a pipe.base.Struct containing the following data.
729 
730  cos_rot_sq : float
731  magnitude of the rotation needed to align the two patterns
732  after their centers are shifted on top of each other.
733  None if rotation test fails.
734  shift_matrix : float array
735  3x3 rotation matrix describing the shift needed to align
736  the source and candidate reference center.
737  None if rotation test fails.
738  """
739 
740  # Make sure the sine is a real number.
741  if cos_shift > 1.0:
742  cos_shift = 1.
743  elif cos_shift < -1.0:
744  cos_shift = -1.
745  sin_shift = np.sqrt(1 - cos_shift ** 2)
746 
747  # If the sine of our shift is zero we only need to use the identity
748  # matrix for the shift. Else we construct the rotation matrix for
749  # shift.
750  if sin_shift > 0:
751  rot_axis = np.cross(src_center, ref_center)
752  rot_axis /= sin_shift
753  shift_matrix = self._create_spherical_rotation_matrix(
754  rot_axis, cos_shift, sin_shift)
755  else:
756  shift_matrix = np.identity(3)
757 
758  # Now that we have our shift we apply it to the src delta vector
759  # and check the rotation.
760  rot_src_delta = np.dot(shift_matrix, src_delta)
761  cos_rot_sq = (np.dot(rot_src_delta, ref_delta) ** 2 /
762  (np.dot(rot_src_delta, rot_src_delta) *
763  np.dot(ref_delta, ref_delta)))
764  # If the rotation isn't in tolerance return None.
765  if cos_rot_sq < max_cos_rot_sq:
766  return pipeBase.Struct(
767  cos_rot_sq=None,
768  shift_matrix=None,)
769  return pipeBase.Struct(
770  cos_rot_sq=cos_rot_sq,
771  shift_matrix=shift_matrix,)
772 
773  def _create_spherical_rotation_matrix(self, rot_axis, cos_rotation,
774  sin_rotion):
775  """Construct a generalized 3D rotation matrix about a given
776  axis.
777 
778  Parameters
779  ----------
780  rot_axis : float array
781  3 vector defining the axis to rotate about.
782  cos_rotation : float
783  cosine of the rotation angle.
784  sin_rotation : float
785  sine of the rotation angle.
786 
787  Return
788  ------
789  float array
790  3x3 spherical, rotation matrix.
791  """
792 
793  rot_cross_matrix = np.array(
794  [[0., -rot_axis[2], rot_axis[1]],
795  [rot_axis[2], 0., -rot_axis[0]],
796  [-rot_axis[1], rot_axis[0], 0.]], dtype=np.float64)
797  shift_matrix = (cos_rotation*np.identity(3) +
798  sin_rotion*rot_cross_matrix +
799  (1. - cos_rotation)*np.outer(rot_axis, rot_axis))
800 
801  return shift_matrix
802 
803  def _create_pattern_spokes(self, src_ctr, src_delta_array, src_dist_array,
804  ref_ctr, ref_ctr_id, ref_delta, ref_dist,
805  ref_delta_array, ref_dist_array,
806  ref_id_array, max_dist_rad, n_match):
807  """ Create the individual spokes that make up the pattern now that the
808  shift and rotation are within tolerance.
809 
810  If we can't create a valid pattern we exit early.
811 
812  Parameters
813  ----------
814  src_ctr : float array
815  3 vector of the source pinwheel center
816  src_delta_array : float array
817  Array of 3 vector deltas between the source center and the pairs
818  that make up the remaining spokes of the pinwheel
819  src_dist_array : float array
820  Array of the distances of each src_delta in the pinwheel
821  ref_ctr : float array
822  3 vector of the candidate reference center
823  ref_ctr_id : int
824  id of the ref_ctr in the master reference array
825  ref_delta : float array
826  3 vector of the first candidate pair of the pinwheel. This is
827  the candidate pair that was matched in the
828  main _construct_pattern_and_shift_rot_matrix loop
829  ref_dist : float
830  Radian distance of the first candidate reference pair
831  ref_delta_array : float array
832  Array of 3 vector deltas that are have the current candidate
833  reference center as part of the pair
834  ref_dist_array : float array
835  Array of vector distances for each of the reference pairs
836  ref_id_array : int array
837  Array of id lookups into the master reference array that our
838  center id object is paired with.
839  max_dist_rad : float
840  Maximum search distance
841  n_match : int
842  Number of source deltas that must be matched into the reference
843  deltas in order to consider this a successful pattern match.
844 
845  Returns
846  -------
847  lsst.pipe.base.Struct
848  The Struct contains the following data:
849 
850  ref_spoke_list : list of ints specifying ids into the master
851  reference array
852  src_spoke_list : list of ints specifying indices into the current
853  source pattern that is being tested.
854  """
855  # Struct where we will be putting our results.
856  output_spokes = pipeBase.Struct(
857  ref_spoke_list=[],
858  src_spoke_list=[],)
859 
860  # Counter for number of spokes we failed to find a reference
861  # candidate for. We break the loop if we haven't found enough.
862  n_fail = 0
863  ref_spoke_list = []
864  src_spoke_list = []
865  # Loop over the source pairs.
866  for src_idx in range(1, len(src_dist_array)):
867  if n_fail > len(src_dist_array) - (n_match - 1):
868  break
869 
870  # Given our length tolerance we can use it to compute a tolerance
871  # on the angle between our spoke.
872  src_sin_tol = (max_dist_rad /
873  (src_dist_array[src_idx] + max_dist_rad))
874  # Test if the small angle approximation will still hold. This is
875  # defined as when sin(theta) ~= theta to within 0.1% of each
876  # other. This also implicitly sets a minimum spoke length that we
877  # can use.
878  if src_sin_tol > 0.0447:
879  n_fail += 1
880  continue
881 
882  # Find the reference pairs that include our candidate pattern
883  # center and sort them in increasing delta
884  ref_dist_idx_array = self._find_candidate_reference_pairs(
885  src_dist_array[src_idx], ref_dist_array, max_dist_rad)
886 
887  # Test the spokes and return the id of the reference object.
888  # Return None if no match is found.
889  ref_id = self._test_spoke(
890  src_ctr, src_delta_array[src_idx], src_dist_array[src_idx],
891  src_delta_array[0], src_dist_array[0], ref_ctr, ref_ctr_id,
892  ref_delta, ref_dist, ref_dist_idx_array, ref_delta_array,
893  ref_dist_array,
894  ref_id_array, src_sin_tol)
895  if ref_id is None:
896  n_fail += 1
897  continue
898 
899  # Append the successful indices to our list. The src_idx needs
900  # an extra iteration to skip the first and second source objects.
901  ref_spoke_list.append(ref_id)
902  src_spoke_list.append(src_idx + 1)
903  # If we found enough reference objects we can return early. This is
904  # n_match - 2 as we already have 2 source objects matched into the
905  # reference data.
906  if len(ref_spoke_list) >= n_match - 2:
907  # Set the struct data and return the struct.
908  output_spokes.ref_spoke_list = ref_spoke_list
909  output_spokes.src_spoke_list = src_spoke_list
910  return output_spokes
911 
912  return output_spokes
913 
914  def _test_spoke(self, src_ctr, src_delta, src_dist, src_ctr_delta,
915  src_ctr_dist, ref_ctr, ref_ctr_id, ref_delta, ref_dist,
916  ref_dist_idx_array, ref_delta_array, ref_dist_array,
917  ref_id_array, src_sin_tol):
918  """Test the opening angle between the first spoke of our pattern
919  for the source object against the reference object.
920 
921  This method makes heavy use of the small angle approximation to perform
922  the comparison.
923 
924  Parameters
925  ----------
926  src_ctr : float array
927  3 vector of the source pinwheel center
928  src_delta : float array
929  3 vector delta from the source center and the source object that
930  makes up the current spoke of the pinwheel we are testing.
931  src_dist : float array
932  Distance of the current spoke we are testing
933  src_ctr_delta : float array
934  3 vector delta between the center of the pattern and the first
935  spoke of the pattern. Used to test compute the opening angle
936  between the current spoke and the first spoke.
937  src_ctr_dist : float
938  Distance between the pairs that make up src_ctr_delta
939  ref_ctr : float array
940  3 vector of the candidate reference center
941  ref_ctr_id : int
942  id lookup of the ref_ctr into the master reference array
943  ref_delta : float array
944  3 vector of the first candidate pair of the pinwheel. That is
945  the candidate pair that was matched in the
946  main _construct_pattern_and_shift_rot_matrix loop
947  ref_dist : float
948  Radian distance of the first candidate reference pair
949  ref_dist_idx_array : int array
950  Indices sorted by the delta distance between the source
951  spoke we are trying to test and the candidate reference
952  spokes.
953  ref_delta_array : float array
954  Array of 3 vector deltas that are have the current candidate
955  reference center as part of the pair
956  ref_dist_array : float array
957  Array of vector distances for each of the reference pairs
958  ref_id_array : int array
959  Array of id lookups into the master reference array that our
960  center id object is paired with.
961  src_sin_tol : float
962  Sine of tolerance allowed between source and reference spoke
963  opening angles.
964 
965  Returns
966  -------
967  None or int
968  If we can not find a candidate spoke we return None else we
969  return an int id into the master reference array.
970  """
971 
972  # Precompute all of the source only cross and dot products so we don't
973  # have to do it for each iteration in the reference loop.
974  cos_theta_src = (np.dot(src_delta, src_ctr_delta) /
975  (src_dist * src_ctr_dist))
976  cross_src = (np.cross(src_delta, src_ctr_delta) /
977  (src_dist * src_ctr_dist))
978  dot_cross_src = np.dot(cross_src, src_ctr)
979 
980  # Loop over our candidate reference objects.
981  for ref_dist_idx in ref_dist_idx_array:
982  # Check the direction of the delta vector.
983  ref_sign = 1
984  if ref_id_array[ref_dist_idx] < ref_ctr_id:
985  ref_sign = -1
986 
987  # Compute the cos between our "center" reference vector and the
988  # current reference candidate.
989  cos_theta_ref = ref_sign * (
990  np.dot(ref_delta_array[ref_dist_idx], ref_delta) /
991  (ref_dist_array[ref_dist_idx] * ref_dist))
992 
993  # Make sure we can safely make the comparison in case
994  # our "center" and candidate vectors are mostly aligned.
995  if cos_theta_ref ** 2 < (1 - src_sin_tol ** 2):
996  cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2 /
997  (1 - cos_theta_ref ** 2))
998  else:
999  cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2 /
1000  src_sin_tol ** 2)
1001  # Test the difference of the cosine of the reference angle against
1002  # the source angle. Assumes that the delta between the two is
1003  # small.
1004  if cos_sq_comparison > src_sin_tol ** 2:
1005  continue
1006 
1007  # The cosine tests the magnitude of the angle but not
1008  # its direction. To do that we need to know the sine as well.
1009  # This cross product calculation does that.
1010  cross_ref = ref_sign * (
1011  np.cross(ref_delta_array[ref_dist_idx], ref_delta) /
1012  (ref_dist_array[ref_dist_idx] * ref_dist))
1013  dot_cross_ref = np.dot(cross_ref, ref_ctr)
1014 
1015  # Check the value of the cos again to make sure that it is not
1016  # near zero.
1017  if abs(cos_theta_src) < src_sin_tol:
1018  sin_comparison = (dot_cross_src - dot_cross_ref) / src_sin_tol
1019  else:
1020  sin_comparison = \
1021  (dot_cross_src - dot_cross_ref) / cos_theta_ref
1022 
1023  if abs(sin_comparison) > src_sin_tol:
1024  continue
1025 
1026  # Return the correct id of the candidate we found.
1027  return ref_id_array[ref_dist_idx]
1028 
1029  return None
1030 
1031  def _create_shift_rot_matrix(self, cos_rot_sq, shift_matrix, src_delta,
1032  ref_ctr, ref_delta):
1033  """ Create the final part of our spherical rotation matrix.
1034 
1035  Parameters
1036  ----------
1037  cos_rot_sq : float
1038  cosine of the rotation needed to align our source and reference
1039  candidate patterns.
1040  shift_matrix : float array
1041  3x3 rotation matrix for shifting the source pattern center on top
1042  of the candidate reference pattern center.
1043  src_delta : float array
1044  3 vector delta of representing the first spoke of the source
1045  pattern
1046  ref_ctr : float array
1047  3 vector on the unit-sphere representing the center of our
1048  reference pattern.
1049  ref_delta : float array
1050  3 vector delta made by the first pair of the reference pattern.
1051 
1052  Returns
1053  -------
1054  lsst.pipe.base.Struct
1055  Struct object containing the following data:
1056 
1057  sin_rot : float sine of the amount of rotation between the
1058  source and reference pattern. We use sine here as it is
1059  signed and tells us the chirality of the rotation.
1060  shift_rot_matrix : float array representing the 3x3 rotation
1061  matrix that takes the source pattern and shifts and rotates
1062  it to align with the reference pattern.
1063  """
1064  cos_rot = np.sqrt(cos_rot_sq)
1065  rot_src_delta = np.dot(shift_matrix, src_delta)
1066  delta_dot_cross = np.dot(np.cross(rot_src_delta, ref_delta), ref_ctr)
1067 
1068  sin_rot = np.sign(delta_dot_cross) * np.sqrt(1 - cos_rot_sq)
1069  rot_matrix = self._create_spherical_rotation_matrix(
1070  ref_ctr, cos_rot, sin_rot)
1071 
1072  shift_rot_matrix = np.dot(rot_matrix, shift_matrix)
1073 
1074  return pipeBase.Struct(
1075  sin_rot=sin_rot,
1076  shift_rot_matrix=shift_rot_matrix,)
1077 
1078  def _intermediate_verify(self, src_pattern, ref_pattern, shift_rot_matrix,
1079  max_dist_rad):
1080  """ Perform an intermediate verify step.
1081 
1082  Rotate the matches references into the source frame and test their
1083  distances against tolerance. Only return true if all points are within
1084  tolerance.
1085 
1086  Parameters
1087  ----------
1088  src_pattern : float array
1089  Array of 3 vectors representing the points that make up our source
1090  pinwheel pattern.
1091  ref_pattern : float array
1092  Array of 3 vectors representing our candidate reference pinwheel
1093  pattern.
1094  shift_rot_matrix : float array
1095  3x3 rotation matrix that takes the source objects and rotates them
1096  onto the frame of the reference objects
1097  max_dist_rad : float
1098  Maximum distance allowed to consider two objects the same.
1099 
1100  Returns
1101  -------
1102  float array or None
1103  Return the fitted shift/rotation matrix if all of the points in our
1104  source pattern are within max_dist_rad of their matched reference
1105  objects. Returns None if this criteria is not satisfied.
1106  """
1107  if len(src_pattern) != len(ref_pattern):
1108  raise ValueError(
1109  "Source pattern length does not match ref pattern.\n"
1110  "\t source pattern len=%i, reference pattern len=%i" %
1111  (len(src_pattern), len(ref_pattern)))
1112 
1114  src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad):
1115  # Now that we know our initial shift and rot matrix is valid we
1116  # want to fit the implied transform using all points from
1117  # our pattern. This is a more robust rotation matrix as our
1118  # initial matrix only used the first 2 points from the source
1119  # pattern to estimate the shift and rotation. The matrix we fit
1120  # are allowed to be non unitary but need to preserve the length of
1121  # the rotated vectors to within the match tolerance.
1122  fit_shift_rot_matrix = least_squares(
1123  _rotation_matrix_chi_sq,
1124  x0=shift_rot_matrix.flatten(),
1125  args=(src_pattern, ref_pattern, max_dist_rad)
1126  ).x.reshape((3, 3))
1127  # Do another verify in case the fits have wondered off.
1129  src_pattern, ref_pattern, fit_shift_rot_matrix,
1130  max_dist_rad):
1131  return fit_shift_rot_matrix
1132 
1133  return None
1134 
1135  def _intermediate_verify_comparison(self, pattern_a, pattern_b,
1136  shift_rot_matrix, max_dist_rad):
1137  """Test the input rotation matrix against one input pattern and
1138  a second one.
1139 
1140  If every point in the pattern after rotation is within a distance of
1141  max_dist_rad to its candidate point in the other pattern, we return
1142  True.
1143 
1144  Parameters
1145  ----------
1146  pattern_a : float array
1147  Array of 3 vectors representing the points that make up our source
1148  pinwheel pattern.
1149  pattern_b : float array
1150  Array of 3 vectors representing our candidate reference pinwheel
1151  pattern.
1152  shift_rot_matrix : float array
1153  3x3 rotation matrix that takes the source objects and rotates them
1154  onto the frame of the reference objects
1155  max_dist_rad : float
1156  Maximum distance allowed to consider two objects the same.
1157 
1158 
1159  Returns
1160  -------
1161  bool
1162  True if all rotated source points are within max_dist_rad of
1163  the candidate references matches.
1164  """
1165  shifted_pattern_a = np.dot(shift_rot_matrix,
1166  pattern_a.transpose()).transpose()
1167  tmp_delta_array = shifted_pattern_a - pattern_b
1168  tmp_dist_array = (tmp_delta_array[:, 0] ** 2 +
1169  tmp_delta_array[:, 1] ** 2 +
1170  tmp_delta_array[:, 2] ** 2)
1171  return np.all(tmp_dist_array < max_dist_rad ** 2)
1172 
1173  def _test_pattern_lengths(self, test_pattern, max_dist_rad):
1174  """ Test that the all vectors in a pattern are unit length within
1175  tolerance.
1176 
1177  This is useful for assuring the non unitary transforms do not contain
1178  too much distortion.
1179  """
1180  dists = (test_pattern[:, 0] ** 2 +
1181  test_pattern[:, 1] ** 2 +
1182  test_pattern[:, 2] ** 2)
1183  return np.all(
1184  np.logical_and((1 - max_dist_rad) ** 2 < dists,
1185  dists < (1 + max_dist_rad) ** 2))
1186 
1187  def _test_rotation_agreement(self, rot_vects, max_dist_rad):
1188  """ Test this rotation against the previous N found and return
1189  the number that a agree within tolerance to where our test
1190  points are.
1191 
1192  Parameters
1193  ----------
1194  rot_vects : float array
1195  Arrays of rotated 3 vectors representing the maximum x, y,
1196  z extent on the unit sphere of the input source objects rotated by
1197  the candidate rotations into the reference frame.
1198  max_dist_rad : float
1199  maximum distance in radians to consider two points "agreeing" on
1200  a rotation
1201 
1202  Returns
1203  -------
1204  int
1205  Number of candidate rotations that agree for all of the rotated
1206  test 3 vectors.
1207  """
1208 
1209  self.log.debug("Comparing pattern %i to previous %i rotations..." %
1210  (rot_vects[-1][-1], len(rot_vects) - 1))
1211 
1212  tot_consent = 0
1213  for rot_idx in range(max((len(rot_vects) - 1), 0)):
1214  tmp_dist_list = []
1215  for vect_idx in range(len(rot_vects[rot_idx]) - 1):
1216  tmp_delta_vect = (rot_vects[rot_idx][vect_idx] -
1217  rot_vects[-1][vect_idx])
1218  tmp_dist_list.append(
1219  np.dot(tmp_delta_vect, tmp_delta_vect))
1220  if np.all(np.array(tmp_dist_list) < max_dist_rad ** 2):
1221  tot_consent += 1
1222  return tot_consent
1223 
1224  def _match_sources(self,
1225  source_array,
1226  shift_rot_matrix):
1227  """ Shift both the reference and source catalog to the the respective
1228  frames and find their nearest neighbor using a kdTree.
1229 
1230  Removes all matches who do not agree when either the reference or
1231  source catalog is rotated. Cuts on a maximum distance are left to an
1232  external function.
1233 
1234  Parameters
1235  ----------
1236  source_array : float array
1237  array of 3 vectors representing the source objects we are trying
1238  to match into the source catalog.
1239  shift_rot_matrix : float array
1240  3x3 rotation matrix that performs the spherical rotation from the
1241  source frame into the reference frame.
1242 
1243  Returns
1244  -------
1245  lsst.pipe.base.Struct
1246  A Struct object containing the following data
1247  matches : a (N, 2) array of integer ids into the source and
1248  reference arrays. Matches are only returned for those that
1249  satisfy the distance and handshake criteria.
1250  distances : float array of the distance between each match in
1251  radians after the shift and rotation is applied.
1252  """
1253  shifted_references = np.dot(
1254  np.linalg.inv(shift_rot_matrix),
1255  self._reference_array.transpose()).transpose()
1256  shifted_sources = np.dot(
1257  shift_rot_matrix,
1258  source_array.transpose()).transpose()
1259 
1260  ref_matches = np.empty((len(shifted_references), 2),
1261  dtype=np.uint32)
1262  src_matches = np.empty((len(shifted_sources), 2),
1263  dtype=np.uint32)
1264 
1265  ref_matches[:, 1] = np.arange(len(shifted_references),
1266  dtype=np.uint32)
1267  src_matches[:, 0] = np.arange(len(shifted_sources),
1268  dtype=np.uint32)
1269 
1270  ref_kdtree = cKDTree(self._reference_array)
1271  src_kdtree = cKDTree(source_array)
1272 
1273  ref_to_src_dist, tmp_ref_to_src_idx = \
1274  src_kdtree.query(shifted_references)
1275  src_to_ref_dist, tmp_src_to_ref_idx = \
1276  ref_kdtree.query(shifted_sources)
1277 
1278  ref_matches[:, 0] = tmp_ref_to_src_idx
1279  src_matches[:, 1] = tmp_src_to_ref_idx
1280 
1281  handshake_mask = self._handshake_match(src_matches, ref_matches)
1282  return pipeBase.Struct(
1283  match_ids=src_matches[handshake_mask],
1284  distances_rad=src_to_ref_dist[handshake_mask],)
1285 
1286  def _handshake_match(self, matches_src, matches_ref):
1287  """Return only those matches where both the source
1288  and reference objects agree they they are each others'
1289  nearest neighbor.
1290  Parameters
1291  ----------
1292  matches_src : int array
1293  (N, 2) int array of nearest neighbor matches between shifted and
1294  rotated reference objects matched into the sources.
1295  matches_ref : int array
1296  (M, 2) int array of nearest neighbor matches between shifted and
1297  rotated source objects matched into the references.
1298  Return
1299  ------
1300  bool array
1301  Return the array positions where the two match catalogs agree.
1302  """
1303  handshake_mask_array = np.zeros(len(matches_src), dtype=np.bool)
1304 
1305  for src_match_idx, match in enumerate(matches_src):
1306  ref_match_idx = np.searchsorted(matches_ref[:, 1], match[1])
1307  if match[0] == matches_ref[ref_match_idx, 0]:
1308  handshake_mask_array[src_match_idx] = True
1309  return handshake_mask_array
def _create_shift_rot_matrix(self, cos_rot_sq, shift_matrix, src_delta, ref_ctr, ref_delta)
def _test_spoke(self, src_ctr, src_delta, src_dist, src_ctr_delta, src_ctr_dist, ref_ctr, ref_ctr_id, ref_delta, ref_dist, ref_dist_idx_array, ref_delta_array, ref_dist_array, ref_id_array, src_sin_tol)
def _construct_pattern_and_shift_rot_matrix(self, src_pattern_array, n_match, max_cos_theta_shift, max_cos_rot_sq, max_dist_rad)
def _test_rotation(self, src_center, ref_center, src_delta, ref_delta, cos_shift, max_cos_rot_sq)
def _intermediate_verify_comparison(self, pattern_a, pattern_b, shift_rot_matrix, max_dist_rad)
def match(self, source_array, n_check, n_match, n_agree, max_n_patterns, max_shift, max_rotation, max_dist, min_matches, pattern_skip_array=None)
def _create_pattern_spokes(self, src_ctr, src_delta_array, src_dist_array, ref_ctr, ref_ctr_id, ref_delta, ref_dist, ref_delta_array, ref_dist_array, ref_id_array, max_dist_rad, n_match)
def _intermediate_verify(self, src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad)