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