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