3 from scipy.optimize
import least_squares
4 from scipy.spatial
import cKDTree
6 import lsst.pipe.base
as pipeBase
9 def _rotation_matrix_chi_sq(flattened_rot_matrix,
13 """Compute the squared differences for least squares fitting. 15 Given a flattened rotation matrix, one N point pattern and another N point 16 pattern to transform into to, compute the squared differences between the 17 points in the two patterns after the rotation. 21 flattened_rot_matrix : `numpy.ndarray`, (9, ) 22 A flattened array representing a 3x3 rotation matrix. The array is 23 flattened to comply with the API of scipy.optimize.least_squares. 24 Flattened elements are [[0, 0], [0, 1], [0, 2], [1, 0]...] 25 pattern_a : `numpy.ndarray`, (N, 3) 26 A array containing N, 3 vectors representing the objects we would like 27 to transform into the frame of pattern_b. 28 pattern_b : `numpy.ndarray`, (N, 3) 29 A array containing N, 3 vectors representing the reference frame we 30 would like to transform pattern_a into. 31 max_dist_rad : `float` 32 The maximum distance allowed from the pattern matching. This value is 33 used as the standard error for the resultant chi values. 37 noralized_diff : `numpy.ndarray`, (9,) 38 Array of differences between the vectors representing of the source 39 pattern rotated into the reference frame and the converse. This is 40 used to minimize in a least squares fitter. 43 rot_matrix = flattened_rot_matrix.reshape((3, 3))
45 rot_pattern_a = np.dot(rot_matrix, pattern_a.transpose()).transpose()
46 diff_pattern_a_to_b = rot_pattern_a - pattern_b
49 return diff_pattern_a_to_b.flatten() / max_dist_rad
53 """Class implementing a pessimistic version of Optimistic Pattern Matcher 54 B (OPMb) from Tabur 2007. See `DMTN-031 <http://ls.st/DMTN-031`_ 58 reference_array : `numpy.ndarray`, (N, 3) 59 spherical points x, y, z of to use as reference objects for 62 Logger for outputting debug info. 66 The class loads and stores the reference object 67 in a convenient data structure for matching any set of source objects that 68 are assumed to contain each other. The pessimistic nature of the algorithm 69 comes from requiring that it discovers at least two patterns that agree on 70 the correct shift and rotation for matching before exiting. The original 71 behavior of OPMb can be recovered simply. Patterns matched between the 72 input datasets are n-spoked pinwheels created from n+1 points. Refer to 73 DMTN #031 for more details. http://github.com/lsst-dm/dmtn-031 83 def _build_distances_and_angles(self):
84 """Create the data structures we will use to search for our pattern 87 Throughout this function and the rest of the class we use id to 88 reference the position in the input reference catalog and index to 89 'index' into the arrays sorted on distance. 103 sub_id_array_list = []
104 sub_dist_array_list = []
112 sub_id_array = np.zeros((self.
_n_reference - 1 - ref_id, 2),
114 sub_id_array[:, 0] = ref_id
115 sub_id_array[:, 1] = np.arange(ref_id + 1, self.
_n_reference,
121 ref_obj).astype(np.float32)
122 sub_dist_array = np.sqrt(sub_delta_array[:, 0] ** 2 +
123 sub_delta_array[:, 1] ** 2 +
124 sub_delta_array[:, 2] ** 2)
128 sub_id_array_list.append(sub_id_array)
129 sub_dist_array_list.append(sub_dist_array)
146 ref_id, sorted_pair_dist_args]
148 ref_id, sorted_pair_dist_args]
151 unsorted_id_array = np.concatenate(sub_id_array_list)
152 unsorted_dist_array = np.concatenate(sub_dist_array_list)
156 sorted_dist_args = unsorted_dist_array.argsort()
157 self.
_dist_array = unsorted_dist_array[sorted_dist_args]
158 self.
_id_array = unsorted_id_array[sorted_dist_args]
162 def match(self, source_array, n_check, n_match, n_agree,
163 max_n_patterns, max_shift, max_rotation, max_dist,
164 min_matches, pattern_skip_array=None):
165 """Match a given source catalog into the loaded reference catalog. 167 Given array of points on the unit sphere and tolerances, we 168 attempt to match a pinwheel like pattern between these input sources 169 and the reference objects this class was created with. This pattern 170 informs of the shift and rotation needed to align the input source 171 objects into the frame of the references. 175 source_array : `numpy.ndarray`, (N, 3) 176 An array of spherical x,y,z coordinates and a magnitude in units 177 of objects having a lower value for sorting. The array should be 180 Number of sources to create a pattern from. Not all objects may be 181 checked if n_match criteria is before looping through all n_check 184 Number of objects to use in constructing a pattern to match. 186 Number of found patterns that must agree on their shift and 187 rotation before exiting. Set this value to 1 to recover the 188 expected behavior of Optimistic Pattern Matcher B. 189 max_n_patters : `int` 190 Number of patterns to create from the input source objects to 191 attempt to match into the reference objects. 193 Maximum allowed shift to match patterns in arcseconds. 194 max_rotation : `float` 195 Maximum allowed rotation between patterns in degrees. 197 Maximum distance in arcseconds allowed between candidate spokes in 198 the source and reference objects. Also sets that maximum distance 199 in the intermediate verify, pattern shift/rotation agreement, and 201 pattern_skip_array : `int` 202 Patterns we would like to skip. This could be due to the pattern 203 being matched on a previous iteration that we now consider invalid. 204 This assumes the ordering of the source objects is the same 205 between different runs of the matcher which, assuming no object 206 has been inserted or the magnitudes have changed, it should be. 210 output_struct : `lsst.pipe.base.Struct` 211 Result struct with components 213 - ``matches`` : (N, 2) array of matched ids for pairs. Empty list if no 214 match found (`numpy.ndarray`, (N, 2) or `list`) 215 - ``distances_rad`` : Radian distances between the matched objects. 216 Empty list if no match found (`numpy.ndarray`, (N,)) 217 - ``pattern_idx``: Index of matched pattern. None if no match found 219 - ``shift`` : Magnitude for the shift between the source and reference 220 objects in arcseconds. None if no match found (`float`). 224 sorted_source_array = source_array[source_array[:, -1].argsort(), :3]
225 n_source = len(sorted_source_array)
228 output_match_struct = pipeBase.Struct(
235 self.
log.warn(
"Source object array is empty. Unable to match. " 249 max_cos_shift = np.cos(np.radians(max_shift / 3600.))
250 max_cos_rot_sq = np.cos(np.radians(max_rotation)) ** 2
251 max_dist_rad = np.radians(max_dist / 3600.)
255 for pattern_idx
in range(np.min((max_n_patterns,
256 n_source - n_match))):
260 if pattern_skip_array
is not None and \
261 np.any(pattern_skip_array == pattern_idx):
263 "Skipping previously matched bad pattern %i..." %
267 pattern = sorted_source_array[
268 pattern_idx: np.min((pattern_idx + n_check, n_source)), :3]
273 construct_return_struct = \
275 pattern, n_match, max_cos_shift, max_cos_rot_sq,
279 if construct_return_struct.ref_candidates
is None or \
280 construct_return_struct.shift_rot_matrix
is None or \
281 construct_return_struct.cos_shift
is None or \
282 construct_return_struct.sin_rot
is None:
286 ref_candidates = construct_return_struct.ref_candidates
287 shift_rot_matrix = construct_return_struct.shift_rot_matrix
288 cos_shift = construct_return_struct.cos_shift
289 sin_rot = construct_return_struct.sin_rot
293 if len(ref_candidates) < n_match:
299 tmp_rot_vect_list = []
300 for test_vect
in test_vectors:
301 tmp_rot_vect_list.append(np.dot(shift_rot_matrix, test_vect))
308 tmp_rot_vect_list.append(pattern_idx)
309 rot_vect_list.append(tmp_rot_vect_list)
321 cut_ids = match_sources_struct.match_ids[
322 match_sources_struct.distances_rad < max_dist_rad]
323 n_matched = len(cut_ids)
326 if n_matched >= min_matches:
332 fit_shift_rot_matrix = least_squares(
333 _rotation_matrix_chi_sq,
334 x0=shift_rot_matrix.flatten(),
335 args=(source_array[cut_ids[:, 0], :3],
340 source_array[:, :3], fit_shift_rot_matrix)
343 shift = np.degrees(np.arccos(cos_shift)) * 3600.
345 self.
log.debug(
"Succeeded after %i patterns." % pattern_idx)
346 self.
log.debug(
"\tShift %.4f arcsec" % shift)
347 self.
log.debug(
"\tRotation: %.4f deg" %
348 np.degrees(np.arcsin(sin_rot)))
351 output_match_struct.match_ids = \
352 match_sources_struct.match_ids
353 output_match_struct.distances_rad = \
354 match_sources_struct.distances_rad
355 output_match_struct.pattern_idx = pattern_idx
356 output_match_struct.shift = shift
357 return output_match_struct
359 self.
log.debug(
"Failed after %i patterns." % (pattern_idx + 1))
360 return output_match_struct
362 def _compute_test_vectors(self, source_array):
363 """Compute spherical 3 vectors at the edges of the x, y, z extent 364 of the input source catalog. 368 source_array : `numpy.ndarray`, (N, 3) 369 array of 3 vectors representing positions on the unit 374 test_vectors : `numpy.ndarray`, (N, 3) 375 Array of vectors representing the maximum extents in x, y, z 376 of the input source array. These are used with the rotations 377 the code finds to test for agreement from different patterns 378 when the code is running in pessimistic mode. 382 if np.any(np.logical_not(np.isfinite(source_array))):
383 self.
log.warn(
"Input source objects contain non-finite values. " 384 "This could end badly.")
385 center_vect = np.nanmean(source_array, axis=0)
389 xbtm_vect = np.array([np.min(source_array[:, 0]), center_vect[1],
390 center_vect[2]], dtype=np.float64)
391 xtop_vect = np.array([np.max(source_array[:, 0]), center_vect[1],
392 center_vect[2]], dtype=np.float64)
393 xbtm_vect /= np.sqrt(np.dot(xbtm_vect, xbtm_vect))
394 xtop_vect /= np.sqrt(np.dot(xtop_vect, xtop_vect))
396 ybtm_vect = np.array([center_vect[0], np.min(source_array[:, 1]),
397 center_vect[2]], dtype=np.float64)
398 ytop_vect = np.array([center_vect[0], np.max(source_array[:, 1]),
399 center_vect[2]], dtype=np.float64)
400 ybtm_vect /= np.sqrt(np.dot(ybtm_vect, ybtm_vect))
401 ytop_vect /= np.sqrt(np.dot(ytop_vect, ytop_vect))
403 zbtm_vect = np.array([center_vect[0], center_vect[1],
404 np.min(source_array[:, 2])], dtype=np.float64)
405 ztop_vect = np.array([center_vect[0], center_vect[1],
406 np.max(source_array[:, 2])], dtype=np.float64)
407 zbtm_vect /= np.sqrt(np.dot(zbtm_vect, zbtm_vect))
408 ztop_vect /= np.sqrt(np.dot(ztop_vect, ztop_vect))
411 return np.array([xbtm_vect, xtop_vect, ybtm_vect, ytop_vect,
412 zbtm_vect, ztop_vect])
414 def _construct_pattern_and_shift_rot_matrix(self, src_pattern_array,
415 n_match, max_cos_theta_shift,
416 max_cos_rot_sq, max_dist_rad):
417 """Test an input source pattern against the reference catalog. 419 Returns the candidate matched patterns and their 420 implied rotation matrices or None. 424 src_pattern_array : `numpy.ndarray`, (N, 3) 425 Sub selection of source 3 vectors to create a pattern from 427 Number of points to attempt to create a pattern from. Must be 428 >= len(src_pattern_array) 429 max_cos_theta_shift : `float` 430 Maximum shift allowed between two patterns' centers. 431 max_cos_rot_sq : `float` 432 Maximum rotation between two patterns that have been shifted 433 to have their centers on top of each other. 434 max_dist_rad : `float` 435 Maximum delta distance allowed between the source and reference 436 pair distances to consider the reference pair a candidate for 437 the source pair. Also sets the tolerance between the opening 438 angles of the spokes when compared to the reference. 442 output_matched_pattern : `lsst.pipe.base.Struct` 443 Result struct with components: 445 - ``ref_candidates`` : ids of the matched pattern in the internal 446 reference_array object (`list` of `int`). 447 - ``src_candidates`` : Pattern ids of the sources matched 449 - ``shift_rot_matrix_src_to_ref`` : 3x3 matrix specifying the full 450 shift and rotation between the reference and source objects. 451 Rotates source into reference frame. `None` if match is not 452 found. (`numpy.ndarray`, (3, 3)) 453 - ``shift_rot_matrix_ref_to_src`` : 3x3 matrix specifying the full 454 shift and rotation of the reference and source objects. Rotates 455 reference into source frame. None if match is not found 456 (`numpy.ndarray`, (3, 3)). 457 - ``cos_shift`` : Magnitude of the shift found between the two 458 patten centers. `None` if match is not found (`float`). 459 - ``sin_rot`` : float value of the rotation to align the already 460 shifted source pattern to the reference pattern. `None` if no match 468 output_matched_pattern = pipeBase.Struct(
471 shift_rot_matrix=
None,
477 src_delta_array = np.empty((len(src_pattern_array) - 1, 3))
478 src_delta_array[:, 0] = (src_pattern_array[1:, 0] -
479 src_pattern_array[0, 0])
480 src_delta_array[:, 1] = (src_pattern_array[1:, 1] -
481 src_pattern_array[0, 1])
482 src_delta_array[:, 2] = (src_pattern_array[1:, 2] -
483 src_pattern_array[0, 2])
484 src_dist_array = np.sqrt(src_delta_array[:, 0]**2 +
485 src_delta_array[:, 1]**2 +
486 src_delta_array[:, 2]**2)
495 for ref_dist_idx
in ref_dist_index_array:
499 tmp_ref_pair_list = self.
_id_array[ref_dist_idx]
500 for pair_idx, ref_id
in enumerate(tmp_ref_pair_list):
501 src_candidates = [0, 1]
503 shift_rot_matrix =
None 510 cos_shift = np.dot(src_pattern_array[0], ref_center)
511 if cos_shift < max_cos_theta_shift:
515 ref_candidates.append(ref_id)
518 ref_candidates.append(
519 tmp_ref_pair_list[1])
523 ref_candidates.append(
524 tmp_ref_pair_list[0])
534 src_pattern_array[0], ref_center, src_delta_array[0],
535 ref_delta, cos_shift, max_cos_rot_sq)
536 if test_rot_struct.cos_rot_sq
is None or \
537 test_rot_struct.proj_ref_ctr_delta
is None or \
538 test_rot_struct.shift_matrix
is None:
542 cos_rot_sq = test_rot_struct.cos_rot_sq
543 proj_ref_ctr_delta = test_rot_struct.proj_ref_ctr_delta
544 shift_matrix = test_rot_struct.shift_matrix
555 src_pattern_array[0], src_delta_array, src_dist_array,
557 tmp_ref_dist_array, tmp_ref_id_array, max_dist_rad,
562 if len(pattern_spoke_struct.ref_spoke_list) < n_match - 2
or \
563 len(pattern_spoke_struct.src_spoke_list) < n_match - 2:
567 ref_candidates.extend(pattern_spoke_struct.ref_spoke_list)
568 src_candidates.extend(pattern_spoke_struct.src_spoke_list)
574 cos_rot_sq, shift_matrix, src_delta_array[0],
578 if shift_rot_struct.sin_rot
is None or \
579 shift_rot_struct.shift_rot_matrix
is None:
583 sin_rot = shift_rot_struct.sin_rot
584 shift_rot_matrix = shift_rot_struct.shift_rot_matrix
593 src_pattern_array[src_candidates],
595 shift_rot_matrix, max_dist_rad)
597 if fit_shift_rot_matrix
is not None:
599 output_matched_pattern.ref_candidates = ref_candidates
600 output_matched_pattern.src_candidates = src_candidates
601 output_matched_pattern.shift_rot_matrix = \
603 output_matched_pattern.cos_shift = cos_shift
604 output_matched_pattern.sin_rot = sin_rot
605 return output_matched_pattern
607 return output_matched_pattern
609 def _find_candidate_reference_pairs(self, src_dist, ref_dist_array,
611 """Wrap numpy.searchsorted to find the range of reference spokes 612 within a spoke distance tolerance of our source spoke. 614 Returns an array sorted from the smallest absolute delta distance 615 between source and reference spoke length. This sorting increases the 616 speed for the pattern search greatly. 621 float value of the distance we would like to search for in 622 the reference array in radians. 623 ref_dist_array : `numpy.ndarray`, (N,) 624 sorted array of distances in radians. 625 max_dist_rad : `float` 626 maximum plus/minus search to find in the reference array in 631 tmp_diff_array : `numpy.ndarray`, (N,) 632 indices lookup into the input ref_dist_array sorted by the 633 difference in value to the src_dist from absolute value 638 start_idx = np.searchsorted(ref_dist_array, src_dist - max_dist_rad)
639 end_idx = np.searchsorted(ref_dist_array, src_dist + max_dist_rad,
643 if start_idx == end_idx:
649 if end_idx > ref_dist_array.shape[0]:
650 end_idx = ref_dist_array.shape[0]
655 tmp_diff_array = np.fabs(ref_dist_array[start_idx:end_idx] - src_dist)
656 return tmp_diff_array.argsort() + start_idx
658 def _test_rotation(self, src_center, ref_center, src_delta, ref_delta,
659 cos_shift, max_cos_rot_sq):
660 """ Test if the rotation implied between the source 661 pattern and reference pattern is within tolerance. To test this 662 we need to create the first part of our spherical rotation matrix 663 which we also return for use later. 667 src_center : `numpy.ndarray`, (N, 3) 669 ref_center : `numpy.ndarray`, (N, 3) 670 3 vector defining the center of the candidate reference pinwheel 672 src_delta : `numpy.ndarray`, (N, 3) 673 3 vector delta between the source pattern center and the end of 675 ref_delta : `numpy.ndarray`, (N, 3) 676 3 vector delta of the candidate matched reference pair 678 Cosine of the angle between the source and reference candidate 680 max_cos_rot_sq : `float` 681 candidate reference pair after shifting the centers on top of each 682 other. The function will return None if the rotation implied is 683 greater than max_cos_rot_sq. 687 result : `lsst.pipe.base.Struct` 688 Result struct with components: 690 - ``cos_rot_sq`` : magnitude of the rotation needed to align the 691 two patterns after their centers are shifted on top of each 692 other. `None` if rotation test fails (`float`). 693 - ``shift_matrix`` : 3x3 rotation matrix describing the shift needed to 694 align the source and candidate reference center. `None` if rotation 695 test fails (`numpy.ndarray`, (N, 3)). 701 elif cos_shift < -1.0:
703 sin_shift = np.sqrt(1 - cos_shift ** 2)
709 rot_axis = np.cross(src_center, ref_center)
710 rot_axis /= sin_shift
712 rot_axis, cos_shift, sin_shift)
714 shift_matrix = np.identity(3)
718 rot_src_delta = np.dot(shift_matrix, src_delta)
719 proj_src_delta = (rot_src_delta -
720 np.dot(rot_src_delta, ref_center) * ref_center)
721 proj_ref_delta = (ref_delta -
722 np.dot(ref_delta, ref_center) * ref_center)
723 cos_rot_sq = (np.dot(proj_src_delta, proj_ref_delta) ** 2 /
724 (np.dot(proj_src_delta, proj_src_delta) *
725 np.dot(proj_ref_delta, proj_ref_delta)))
727 if cos_rot_sq < max_cos_rot_sq:
728 return pipeBase.Struct(
730 proj_ref_ctr_delta=
None,
734 return pipeBase.Struct(
735 cos_rot_sq=cos_rot_sq,
736 proj_ref_ctr_delta=proj_ref_delta,
737 shift_matrix=shift_matrix,)
739 def _create_spherical_rotation_matrix(self, rot_axis, cos_rotation,
741 """Construct a generalized 3D rotation matrix about a given 746 rot_axis : `numpy.ndarray`, (3,) 747 3 vector defining the axis to rotate about. 748 cos_rotation : `float` 749 cosine of the rotation angle. 750 sin_rotation : `float` 751 sine of the rotation angle. 755 shift_matrix : `numpy.ndarray`, (3, 3) 756 3x3 spherical, rotation matrix. 759 rot_cross_matrix = np.array(
760 [[0., -rot_axis[2], rot_axis[1]],
761 [rot_axis[2], 0., -rot_axis[0]],
762 [-rot_axis[1], rot_axis[0], 0.]], dtype=np.float64)
763 shift_matrix = (cos_rotation*np.identity(3) +
764 sin_rotion*rot_cross_matrix +
765 (1. - cos_rotation)*np.outer(rot_axis, rot_axis))
769 def _create_pattern_spokes(self, src_ctr, src_delta_array, src_dist_array,
770 ref_ctr, ref_ctr_id, proj_ref_ctr_delta,
771 ref_dist_array, ref_id_array, max_dist_rad,
773 """ Create the individual spokes that make up the pattern now that the 774 shift and rotation are within tolerance. 776 If we can't create a valid pattern we exit early. 780 src_ctr : `numpy.ndarray`, (3,) 781 3 vector of the source pinwheel center 782 src_delta_array : `numpy.ndarray`, (N, 3) 783 Array of 3 vector deltas between the source center and the pairs 784 that make up the remaining spokes of the pinwheel 785 src_dist_array : `numpy.ndarray`, (N, 3) 786 Array of the distances of each src_delta in the pinwheel 787 ref_ctr : `numpy.ndarray`, (3,) 788 3 vector of the candidate reference center 790 id of the ref_ctr in the master reference array 791 proj_ref_ctr_delta : `numpy.ndarray`, (3,) 792 Plane projected 3 vector formed from the center point of the 793 candidate pin-wheel and the second point in the pattern to create 794 the first spoke pair. This is the candidate pair that was matched 795 in the main _construct_pattern_and_shift_rot_matrix loop 796 ref_dist_array : `numpy.ndarray`, (N,) 797 Array of vector distances for each of the reference pairs 798 ref_id_array : `numpy.ndarray`, (N,) 799 Array of id lookups into the master reference array that our 800 center id object is paired with. 801 max_dist_rad : `float` 802 Maximum search distance 804 Number of source deltas that must be matched into the reference 805 deltas in order to consider this a successful pattern match. 809 output_spokes : `lsst.pipe.base.Struct` 810 Result struct with components: 812 - ``ref_spoke_list`` : list of ints specifying ids into the master 813 reference array (`list` of `int`). 814 - ``src_spoke_list`` : list of ints specifying indices into the 815 current source pattern that is being tested (`list` of `int`). 818 output_spokes = pipeBase.Struct(
830 proj_src_ctr_delta = (src_delta_array[0] -
831 np.dot(src_delta_array[0], src_ctr) * src_ctr)
832 proj_src_ctr_dist_sq = np.dot(proj_src_ctr_delta, proj_src_ctr_delta)
835 proj_ref_ctr_dist_sq = np.dot(proj_ref_ctr_delta, proj_ref_ctr_delta)
838 for src_idx
in range(1, len(src_dist_array)):
839 if n_fail > len(src_dist_array) - (n_match - 1):
844 src_sin_tol = (max_dist_rad /
845 (src_dist_array[src_idx] + max_dist_rad))
852 if src_sin_tol > max_sin_tol:
853 src_sin_tol = max_sin_tol
858 src_delta_array[src_idx] -
859 np.dot(src_delta_array[src_idx], src_ctr) * src_ctr)
860 geom_dist_src = np.sqrt(
861 np.dot(proj_src_delta, proj_src_delta) *
862 proj_src_ctr_dist_sq)
865 cos_theta_src = (np.dot(proj_src_delta, proj_src_ctr_delta) /
867 cross_src = (np.cross(proj_src_delta, proj_src_ctr_delta) /
869 sin_theta_src = np.dot(cross_src, src_ctr)
874 src_dist_array[src_idx], ref_dist_array, max_dist_rad)
884 proj_ref_ctr_dist_sq,
894 ref_spoke_list.append(ref_id)
895 src_spoke_list.append(src_idx + 1)
899 if len(ref_spoke_list) >= n_match - 2:
901 output_spokes.ref_spoke_list = ref_spoke_list
902 output_spokes.src_spoke_list = src_spoke_list
907 def _test_spoke(self, cos_theta_src, sin_theta_src, ref_ctr, ref_ctr_id,
908 proj_ref_ctr_delta, proj_ref_ctr_dist_sq,
909 ref_dist_idx_array, ref_id_array, src_sin_tol):
910 """Test the opening angle between the first spoke of our pattern 911 for the source object against the reference object. 913 This method makes heavy use of the small angle approximation to perform 918 cos_theta_src : `float` 919 Cosine of the angle between the current candidate source spoke and 921 sin_theta_src : `float` 922 Sine of the angle between the current candidate source spoke and 924 ref_ctr : `numpy.ndarray`, (3,) 925 3 vector of the candidate reference center 927 id lookup of the ref_ctr into the master reference array 928 proj_ref_ctr_delta : `float` 929 Plane projected first spoke in the reference pattern using the 930 pattern center as normal. 931 proj_ref_ctr_dist_sq : `float` 932 Squared length of the projected vector. 933 ref_dist_idx_array : `numpy.ndarray`, (N,) 934 Indices sorted by the delta distance between the source 935 spoke we are trying to test and the candidate reference 937 ref_id_array : `numpy.ndarray`, (N,) 938 Array of id lookups into the master reference array that our 939 center id object is paired with. 940 src_sin_tol : `float` 941 Sine of tolerance allowed between source and reference spoke 947 If we can not find a candidate spoke we return `None` else we 948 return an int id into the master reference array. 952 for ref_dist_idx
in ref_dist_idx_array:
958 proj_ref_delta = ref_delta - np.dot(ref_delta, ref_ctr) * ref_ctr
959 geom_dist_ref = np.sqrt(proj_ref_ctr_dist_sq *
960 np.dot(proj_ref_delta, proj_ref_delta))
961 cos_theta_ref = (np.dot(proj_ref_delta, proj_ref_ctr_delta) /
966 if cos_theta_ref ** 2 < (1 - src_sin_tol ** 2):
967 cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2 /
968 (1 - cos_theta_ref ** 2))
970 cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2 /
975 if cos_sq_comparison > src_sin_tol ** 2:
981 cross_ref = (np.cross(proj_ref_delta, proj_ref_ctr_delta) /
983 sin_theta_ref = np.dot(cross_ref, ref_ctr)
987 if abs(cos_theta_src) < src_sin_tol:
988 sin_comparison = (sin_theta_src - sin_theta_ref) / src_sin_tol
991 (sin_theta_src - sin_theta_ref) / cos_theta_ref
993 if abs(sin_comparison) > src_sin_tol:
997 return ref_id_array[ref_dist_idx]
1001 def _create_shift_rot_matrix(self, cos_rot_sq, shift_matrix, src_delta,
1002 ref_ctr, ref_delta):
1003 """ Create the final part of our spherical rotation matrix. 1007 cos_rot_sq : `float` 1008 cosine of the rotation needed to align our source and reference 1010 shift_matrix : `numpy.ndarray`, (3, 3) 1011 3x3 rotation matrix for shifting the source pattern center on top 1012 of the candidate reference pattern center. 1013 src_delta : `numpy.ndarray`, (3,) 1014 3 vector delta of representing the first spoke of the source 1016 ref_ctr : `numpy.ndarray`, (3,) 1017 3 vector on the unit-sphere representing the center of our 1019 ref_delta : `numpy.ndarray`, (3,) 1020 3 vector delta made by the first pair of the reference pattern. 1024 result : `lsst.pipe.base.Struct` 1025 Result struct with components: 1027 - ``sin_rot`` : float sine of the amount of rotation between the 1028 source and reference pattern. We use sine here as it is 1029 signed and tells us the chirality of the rotation (`float`). 1030 - ``shift_rot_matrix`` : float array representing the 3x3 rotation 1031 matrix that takes the source pattern and shifts and rotates 1032 it to align with the reference pattern (`numpy.ndarray`, (3,3)). 1034 cos_rot = np.sqrt(cos_rot_sq)
1035 rot_src_delta = np.dot(shift_matrix, src_delta)
1036 delta_dot_cross = np.dot(np.cross(rot_src_delta, ref_delta), ref_ctr)
1038 sin_rot = np.sign(delta_dot_cross) * np.sqrt(1 - cos_rot_sq)
1040 ref_ctr, cos_rot, sin_rot)
1042 shift_rot_matrix = np.dot(rot_matrix, shift_matrix)
1044 return pipeBase.Struct(
1046 shift_rot_matrix=shift_rot_matrix,)
1048 def _intermediate_verify(self, src_pattern, ref_pattern, shift_rot_matrix,
1050 """ Perform an intermediate verify step. 1052 Rotate the matches references into the source frame and test their 1053 distances against tolerance. Only return true if all points are within 1058 src_pattern : `numpy.ndarray`, (N,3) 1059 Array of 3 vectors representing the points that make up our source 1061 ref_pattern : `numpy.ndarray`, (N,3) 1062 Array of 3 vectors representing our candidate reference pinwheel 1064 shift_rot_matrix : `numpy.ndarray`, (3,3) 1065 3x3 rotation matrix that takes the source objects and rotates them 1066 onto the frame of the reference objects 1067 max_dist_rad : `float` 1068 Maximum distance allowed to consider two objects the same. 1072 fit_shift_rot_matrix : `numpy.ndarray`, (3,3) 1073 Return the fitted shift/rotation matrix if all of the points in our 1074 source pattern are within max_dist_rad of their matched reference 1075 objects. Returns None if this criteria is not satisfied. 1077 if len(src_pattern) != len(ref_pattern):
1079 "Source pattern length does not match ref pattern.\n" 1080 "\t source pattern len=%i, reference pattern len=%i" %
1081 (len(src_pattern), len(ref_pattern)))
1084 src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad):
1092 fit_shift_rot_matrix = least_squares(
1093 _rotation_matrix_chi_sq,
1094 x0=shift_rot_matrix.flatten(),
1095 args=(src_pattern, ref_pattern, max_dist_rad)
1099 src_pattern, ref_pattern, fit_shift_rot_matrix,
1101 return fit_shift_rot_matrix
1105 def _intermediate_verify_comparison(self, pattern_a, pattern_b,
1106 shift_rot_matrix, max_dist_rad):
1107 """Test the input rotation matrix against one input pattern and 1110 If every point in the pattern after rotation is within a distance of 1111 max_dist_rad to its candidate point in the other pattern, we return 1116 pattern_a : `numpy.ndarray`, (N,3) 1117 Array of 3 vectors representing the points that make up our source 1119 pattern_b : `numpy.ndarray`, (N,3) 1120 Array of 3 vectors representing our candidate reference pinwheel 1122 shift_rot_matrix : `numpy.ndarray`, (3,3) 1123 3x3 rotation matrix that takes the source objects and rotates them 1124 onto the frame of the reference objects 1125 max_dist_rad : `float` 1126 Maximum distance allowed to consider two objects the same. 1132 True if all rotated source points are within max_dist_rad of 1133 the candidate references matches. 1135 shifted_pattern_a = np.dot(shift_rot_matrix,
1136 pattern_a.transpose()).transpose()
1137 tmp_delta_array = shifted_pattern_a - pattern_b
1138 tmp_dist_array = (tmp_delta_array[:, 0] ** 2 +
1139 tmp_delta_array[:, 1] ** 2 +
1140 tmp_delta_array[:, 2] ** 2)
1141 return np.all(tmp_dist_array < max_dist_rad ** 2)
1143 def _test_pattern_lengths(self, test_pattern, max_dist_rad):
1144 """ Test that the all vectors in a pattern are unit length within 1147 This is useful for assuring the non unitary transforms do not contain 1148 too much distortion. 1152 test_pattern : `numpy.ndarray`, (N, 3) 1153 Test vectors at the maximum and minimum x, y, z extents. 1154 max_dist_rad : `float` 1155 maximum distance in radians to consider two points "agreeing" on 1163 dists = (test_pattern[:, 0] ** 2 +
1164 test_pattern[:, 1] ** 2 +
1165 test_pattern[:, 2] ** 2)
1167 np.logical_and((1 - max_dist_rad) ** 2 < dists,
1168 dists < (1 + max_dist_rad) ** 2))
1170 def _test_rotation_agreement(self, rot_vects, max_dist_rad):
1171 """ Test this rotation against the previous N found and return 1172 the number that a agree within tolerance to where our test 1177 rot_vects : `numpy.ndarray`, (N, 3) 1178 Arrays of rotated 3 vectors representing the maximum x, y, 1179 z extent on the unit sphere of the input source objects rotated by 1180 the candidate rotations into the reference frame. 1181 max_dist_rad : `float` 1182 maximum distance in radians to consider two points "agreeing" on 1188 Number of candidate rotations that agree for all of the rotated 1192 self.
log.debug(
"Comparing pattern %i to previous %i rotations..." %
1193 (rot_vects[-1][-1], len(rot_vects) - 1))
1196 for rot_idx
in range(max((len(rot_vects) - 1), 0)):
1198 for vect_idx
in range(len(rot_vects[rot_idx]) - 1):
1199 tmp_delta_vect = (rot_vects[rot_idx][vect_idx] -
1200 rot_vects[-1][vect_idx])
1201 tmp_dist_list.append(
1202 np.dot(tmp_delta_vect, tmp_delta_vect))
1203 if np.all(np.array(tmp_dist_list) < max_dist_rad ** 2):
1207 def _match_sources(self,
1210 """ Shift both the reference and source catalog to the the respective 1211 frames and find their nearest neighbor using a kdTree. 1213 Removes all matches who do not agree when either the reference or 1214 source catalog is rotated. Cuts on a maximum distance are left to an 1219 source_array : `numpy.ndarray`, (N, 3) 1220 array of 3 vectors representing the source objects we are trying 1221 to match into the source catalog. 1222 shift_rot_matrix : `numpy.ndarray`, (3, 3) 1223 3x3 rotation matrix that performs the spherical rotation from the 1224 source frame into the reference frame. 1228 results : `lsst.pipe.base.Struct` 1229 Result struct with components: 1231 - ``matches`` : array of integer ids into the source and 1232 reference arrays. Matches are only returned for those that 1233 satisfy the distance and handshake criteria 1234 (`numpy.ndarray`, (N, 2)). 1235 - ``distances`` : Distances between each match in radians after 1236 the shift and rotation is applied (`numpy.ndarray`, (N)). 1238 shifted_references = np.dot(
1239 np.linalg.inv(shift_rot_matrix),
1241 shifted_sources = np.dot(
1243 source_array.transpose()).transpose()
1245 ref_matches = np.empty((len(shifted_references), 2),
1247 src_matches = np.empty((len(shifted_sources), 2),
1250 ref_matches[:, 1] = np.arange(len(shifted_references),
1252 src_matches[:, 0] = np.arange(len(shifted_sources),
1256 src_kdtree = cKDTree(source_array)
1258 ref_to_src_dist, tmp_ref_to_src_idx = \
1259 src_kdtree.query(shifted_references)
1260 src_to_ref_dist, tmp_src_to_ref_idx = \
1261 ref_kdtree.query(shifted_sources)
1263 ref_matches[:, 0] = tmp_ref_to_src_idx
1264 src_matches[:, 1] = tmp_src_to_ref_idx
1267 return pipeBase.Struct(
1268 match_ids=src_matches[handshake_mask],
1269 distances_rad=src_to_ref_dist[handshake_mask],)
1271 def _handshake_match(self, matches_src, matches_ref):
1272 """Return only those matches where both the source 1273 and reference objects agree they they are each others' 1278 matches_src : `numpy.ndarray`, (N, 2) 1279 int array of nearest neighbor matches between shifted and 1280 rotated reference objects matched into the sources. 1281 matches_ref : `numpy.ndarray`, (N, 2) 1282 int array of nearest neighbor matches between shifted and 1283 rotated source objects matched into the references. 1286 handshake_mask_array : `numpy.ndarray`, (N,) 1287 Return the array positions where the two match catalogs agree. 1289 handshake_mask_array = np.zeros(len(matches_src), dtype=np.bool)
1291 for src_match_idx, match
in enumerate(matches_src):
1292 ref_match_idx = np.searchsorted(matches_ref[:, 1], match[1])
1293 if match[0] == matches_ref[ref_match_idx, 0]:
1294 handshake_mask_array[src_match_idx] =
True 1295 return handshake_mask_array
def _create_shift_rot_matrix(self, cos_rot_sq, shift_matrix, src_delta, ref_ctr, ref_delta)
def _test_pattern_lengths(self, test_pattern, max_dist_rad)
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 _match_sources(self, source_array, shift_rot_matrix)
def __init__(self, reference_array, log)
def _handshake_match(self, matches_src, matches_ref)
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 _build_distances_and_angles(self)
def _test_rotation_agreement(self, rot_vects, max_dist_rad)
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 _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 _compute_test_vectors(self, source_array)
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 _find_candidate_reference_pairs(self, src_dist, ref_dist_array, max_dist_rad)
def _create_spherical_rotation_matrix(self, rot_axis, cos_rotation, sin_rotion)
def _intermediate_verify(self, src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad)