Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1 

2import numpy as np 

3from scipy.optimize import least_squares 

4from scipy.spatial import cKDTree 

5from scipy.stats import sigmaclip 

6 

7import lsst.pipe.base as pipeBase 

8 

9 

10def _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 

53class PessimisticPatternMatcherB: 

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 

79 self._n_reference = len(self._reference_array) 

80 self.log = log 

81 

82 self._build_distances_and_angles() 

83 

84 def _build_distances_and_angles(self): 

85 """Create the data structures we will use to search for our pattern 

86 match in. 

87 

88 Throughout this function and the rest of the class we use id to 

89 reference the position in the input reference catalog and index to 

90 'index' into the arrays sorted on distance. 

91 """ 

92 

93 # Initialize the arrays we will need for quick look up of pairs once 

94 # have a candidate spoke center. 

95 self._pair_id_array = np.empty( 

96 (self._n_reference, self._n_reference - 1), 

97 dtype=np.uint16) 

98 self._pair_dist_array = np.empty( 

99 (self._n_reference, self._n_reference - 1), 

100 dtype=np.float32) 

101 

102 # Create empty lists to temporarily store our pair information per 

103 # reference object. These will be concatenated into our final arrays. 

104 sub_id_array_list = [] 

105 sub_dist_array_list = [] 

106 

107 # Loop over reference objects storing pair distances and ids. 

108 for ref_id, ref_obj in enumerate(self._reference_array): 

109 

110 # Reserve and fill the ids of each reference object pair. 

111 # 16 bit is safe for the id array as the catalog input from 

112 # MatchPessimisticB is limited to a max length of 2 ** 16. 

113 sub_id_array = np.zeros((self._n_reference - 1 - ref_id, 2), 

114 dtype=np.uint16) 

115 sub_id_array[:, 0] = ref_id 

116 sub_id_array[:, 1] = np.arange(ref_id + 1, self._n_reference, 

117 dtype=np.uint16) 

118 

119 # Compute the vector deltas for each pair of reference objects. 

120 # Compute and store the distances. 

121 sub_delta_array = (self._reference_array[ref_id + 1:, :] 

122 - ref_obj).astype(np.float32) 

123 sub_dist_array = np.sqrt(sub_delta_array[:, 0] ** 2 

124 + sub_delta_array[:, 1] ** 2 

125 + sub_delta_array[:, 2] ** 2) 

126 

127 # Append to our arrays to the output lists for later 

128 # concatenation. 

129 sub_id_array_list.append(sub_id_array) 

130 sub_dist_array_list.append(sub_dist_array) 

131 

132 # Fill the pair look up arrays row wise and then column wise. 

133 self._pair_id_array[ref_id, ref_id:] = sub_id_array[:, 1] 

134 self._pair_dist_array[ref_id, ref_id:] = sub_dist_array 

135 

136 # Don't fill the array column wise if we are on the last as 

137 # these pairs have already been filled by previous 

138 # iterations. 

139 if ref_id < self._n_reference - 1: 

140 self._pair_id_array[ref_id + 1:, ref_id] = ref_id 

141 self._pair_dist_array[ref_id + 1:, ref_id] = sub_dist_array 

142 

143 # Sort each row on distance for fast look up of pairs given 

144 # the id of one of the objects in the pair. 

145 sorted_pair_dist_args = self._pair_dist_array[ref_id, :].argsort() 

146 self._pair_dist_array[ref_id, :] = self._pair_dist_array[ 

147 ref_id, sorted_pair_dist_args] 

148 self._pair_id_array[ref_id, :] = self._pair_id_array[ 

149 ref_id, sorted_pair_dist_args] 

150 

151 # Concatenate our arrays together. 

152 unsorted_id_array = np.concatenate(sub_id_array_list) 

153 unsorted_dist_array = np.concatenate(sub_dist_array_list) 

154 

155 # Sort each array on the pair distances for the initial 

156 # optimistic pattern matcher lookup. 

157 sorted_dist_args = unsorted_dist_array.argsort() 

158 self._dist_array = unsorted_dist_array[sorted_dist_args] 

159 self._id_array = unsorted_id_array[sorted_dist_args] 

160 

161 return None 

162 

163 def match(self, source_array, n_check, n_match, n_agree, 

164 max_n_patterns, max_shift, max_rotation, max_dist, 

165 min_matches, pattern_skip_array=None): 

166 """Match a given source catalog into the loaded reference catalog. 

167 

168 Given array of points on the unit sphere and tolerances, we 

169 attempt to match a pinwheel like pattern between these input sources 

170 and the reference objects this class was created with. This pattern 

171 informs of the shift and rotation needed to align the input source 

172 objects into the frame of the references. 

173 

174 Parameters 

175 ---------- 

176 source_array : `numpy.ndarray`, (N, 3) 

177 An array of spherical x,y,z coordinates and a magnitude in units 

178 of objects having a lower value for sorting. The array should be 

179 of shape (N, 4). 

180 n_check : `int` 

181 Number of sources to create a pattern from. Not all objects may be 

182 checked if n_match criteria is before looping through all n_check 

183 objects. 

184 n_match : `int` 

185 Number of objects to use in constructing a pattern to match. 

186 n_agree : `int` 

187 Number of found patterns that must agree on their shift and 

188 rotation before exiting. Set this value to 1 to recover the 

189 expected behavior of Optimistic Pattern Matcher B. 

190 max_n_patters : `int` 

191 Number of patterns to create from the input source objects to 

192 attempt to match into the reference objects. 

193 max_shift : `float` 

194 Maximum allowed shift to match patterns in arcseconds. 

195 max_rotation : `float` 

196 Maximum allowed rotation between patterns in degrees. 

197 max_dist : `float` 

198 Maximum distance in arcseconds allowed between candidate spokes in 

199 the source and reference objects. Also sets that maximum distance 

200 in the intermediate verify, pattern shift/rotation agreement, and 

201 final verify steps. 

202 pattern_skip_array : `int` 

203 Patterns we would like to skip. This could be due to the pattern 

204 being matched on a previous iteration that we now consider invalid. 

205 This assumes the ordering of the source objects is the same 

206 between different runs of the matcher which, assuming no object 

207 has been inserted or the magnitudes have changed, it should be. 

208 

209 Returns 

210 ------- 

211 output_struct : `lsst.pipe.base.Struct` 

212 Result struct with components 

213 

214 - ``matches`` : (N, 2) array of matched ids for pairs. Empty list if no 

215 match found (`numpy.ndarray`, (N, 2) or `list`) 

216 - ``distances_rad`` : Radian distances between the matched objects. 

217 Empty list if no match found (`numpy.ndarray`, (N,)) 

218 - ``pattern_idx``: Index of matched pattern. None if no match found 

219 (`int`). 

220 - ``shift`` : Magnitude for the shift between the source and reference 

221 objects in arcseconds. None if no match found (`float`). 

222 """ 

223 

224 # Given our input source_array we sort on magnitude. 

225 sorted_source_array = source_array[source_array[:, -1].argsort(), :3] 

226 n_source = len(sorted_source_array) 

227 

228 # Initialize output struct. 

229 output_match_struct = pipeBase.Struct( 

230 match_ids=[], 

231 distances_rad=[], 

232 pattern_idx=None, 

233 shift=None, 

234 max_dist_rad=None,) 

235 

236 if n_source <= 0: 

237 self.log.warn("Source object array is empty. Unable to match. " 

238 "Exiting matcher.") 

239 return None 

240 

241 # To test if the shifts and rotations we find agree with each other, 

242 # we first create two test points situated at the top and bottom of 

243 # where the z axis on the sphere bisects the source catalog. 

244 test_vectors = self._compute_test_vectors(source_array[:, :3]) 

245 

246 # We now create an empty list of our resultant rotated vectors to 

247 # compare the different rotations we find. 

248 rot_vect_list = [] 

249 

250 # Convert the tolerances to values we will use in the code. 

251 max_cos_shift = np.cos(np.radians(max_shift / 3600.)) 

252 max_cos_rot_sq = np.cos(np.radians(max_rotation)) ** 2 

253 max_dist_rad = np.radians(max_dist / 3600.) 

254 

255 # Loop through the sources from brightest to faintest, grabbing a 

256 # chunk of n_check each time. 

257 for pattern_idx in range(np.min((max_n_patterns, 

258 n_source - n_match))): 

259 

260 # If this pattern is one that we matched on the past but we 

261 # now want to skip, we do so here. 

262 if pattern_skip_array is not None and \ 

263 np.any(pattern_skip_array == pattern_idx): 

264 self.log.debug( 

265 "Skipping previously matched bad pattern %i..." % 

266 pattern_idx) 

267 continue 

268 # Grab the sources to attempt to create this pattern. 

269 pattern = sorted_source_array[ 

270 pattern_idx: np.min((pattern_idx + n_check, n_source)), :3] 

271 

272 # Construct a pattern given the number of points defining the 

273 # pattern complexity. This is the start of the primary tests to 

274 # match our source pattern into the reference objects. 

275 construct_return_struct = \ 

276 self._construct_pattern_and_shift_rot_matrix( 

277 pattern, n_match, max_cos_shift, max_cos_rot_sq, 

278 max_dist_rad) 

279 

280 # Our struct is None if we could not match the pattern. 

281 if construct_return_struct.ref_candidates is None or \ 

282 construct_return_struct.shift_rot_matrix is None or \ 

283 construct_return_struct.cos_shift is None or \ 

284 construct_return_struct.sin_rot is None: 

285 continue 

286 

287 # Grab the output data from the Struct object. 

288 ref_candidates = construct_return_struct.ref_candidates 

289 shift_rot_matrix = construct_return_struct.shift_rot_matrix 

290 cos_shift = construct_return_struct.cos_shift 

291 sin_rot = construct_return_struct.sin_rot 

292 

293 # If we didn't match enough candidates we continue to the next 

294 # pattern. 

295 if len(ref_candidates) < n_match: 

296 continue 

297 

298 # Now that we know our pattern and shift/rotation are valid we 

299 # store the the rotated versions of our test points for later 

300 # use. 

301 tmp_rot_vect_list = [] 

302 for test_vect in test_vectors: 

303 tmp_rot_vect_list.append(np.dot(shift_rot_matrix, test_vect)) 

304 # Since our test point are in the source frame, we can test if 

305 # their lengths are mostly preserved under the transform. 

306 if not self._test_pattern_lengths(np.array(tmp_rot_vect_list), 

307 max_dist_rad): 

308 continue 

309 

310 tmp_rot_vect_list.append(pattern_idx) 

311 rot_vect_list.append(tmp_rot_vect_list) 

312 

313 # Test if we have enough rotations, which agree, or if we 

314 # are in optimistic mode. 

315 if self._test_rotation_agreement(rot_vect_list, max_dist_rad) < \ 

316 n_agree - 1: 

317 continue 

318 

319 # Run the final verify step. 

320 match_struct = self._final_verify(source_array[:, :3], 

321 shift_rot_matrix, 

322 max_dist_rad, 

323 min_matches) 

324 if match_struct.match_ids is None or \ 

325 match_struct.distances_rad is None or \ 

326 match_struct.max_dist_rad is None: 

327 continue 

328 

329 # Convert the observed shift to arcseconds 

330 shift = np.degrees(np.arccos(cos_shift)) * 3600. 

331 # Print information to the logger. 

332 self.log.debug("Succeeded after %i patterns." % pattern_idx) 

333 self.log.debug("\tShift %.4f arcsec" % shift) 

334 self.log.debug("\tRotation: %.4f deg" % 

335 np.degrees(np.arcsin(sin_rot))) 

336 

337 # Fill the struct and return. 

338 output_match_struct.match_ids = \ 

339 match_struct.match_ids 

340 output_match_struct.distances_rad = \ 

341 match_struct.distances_rad 

342 output_match_struct.pattern_idx = pattern_idx 

343 output_match_struct.shift = shift 

344 output_match_struct.max_dist_rad = match_struct.max_dist_rad 

345 return output_match_struct 

346 

347 self.log.debug("Failed after %i patterns." % (pattern_idx + 1)) 

348 return output_match_struct 

349 

350 def _compute_test_vectors(self, source_array): 

351 """Compute spherical 3 vectors at the edges of the x, y, z extent 

352 of the input source catalog. 

353 

354 Parameters 

355 ---------- 

356 source_array : `numpy.ndarray`, (N, 3) 

357 array of 3 vectors representing positions on the unit 

358 sphere. 

359 

360 Returns 

361 ------- 

362 test_vectors : `numpy.ndarray`, (N, 3) 

363 Array of vectors representing the maximum extents in x, y, z 

364 of the input source array. These are used with the rotations 

365 the code finds to test for agreement from different patterns 

366 when the code is running in pessimistic mode. 

367 """ 

368 

369 # Get the center of source_array. 

370 if np.any(np.logical_not(np.isfinite(source_array))): 

371 self.log.warn("Input source objects contain non-finite values. " 

372 "This could end badly.") 

373 center_vect = np.nanmean(source_array, axis=0) 

374 

375 # So that our rotation test works over the full sky we compute 

376 # the max extent in each Cartesian direction x,y,z. 

377 xbtm_vect = np.array([np.min(source_array[:, 0]), center_vect[1], 

378 center_vect[2]], dtype=np.float64) 

379 xtop_vect = np.array([np.max(source_array[:, 0]), center_vect[1], 

380 center_vect[2]], dtype=np.float64) 

381 xbtm_vect /= np.sqrt(np.dot(xbtm_vect, xbtm_vect)) 

382 xtop_vect /= np.sqrt(np.dot(xtop_vect, xtop_vect)) 

383 

384 ybtm_vect = np.array([center_vect[0], np.min(source_array[:, 1]), 

385 center_vect[2]], dtype=np.float64) 

386 ytop_vect = np.array([center_vect[0], np.max(source_array[:, 1]), 

387 center_vect[2]], dtype=np.float64) 

388 ybtm_vect /= np.sqrt(np.dot(ybtm_vect, ybtm_vect)) 

389 ytop_vect /= np.sqrt(np.dot(ytop_vect, ytop_vect)) 

390 

391 zbtm_vect = np.array([center_vect[0], center_vect[1], 

392 np.min(source_array[:, 2])], dtype=np.float64) 

393 ztop_vect = np.array([center_vect[0], center_vect[1], 

394 np.max(source_array[:, 2])], dtype=np.float64) 

395 zbtm_vect /= np.sqrt(np.dot(zbtm_vect, zbtm_vect)) 

396 ztop_vect /= np.sqrt(np.dot(ztop_vect, ztop_vect)) 

397 

398 # Return our list of vectors for later rotation testing. 

399 return np.array([xbtm_vect, xtop_vect, ybtm_vect, ytop_vect, 

400 zbtm_vect, ztop_vect]) 

401 

402 def _construct_pattern_and_shift_rot_matrix(self, src_pattern_array, 

403 n_match, max_cos_theta_shift, 

404 max_cos_rot_sq, max_dist_rad): 

405 """Test an input source pattern against the reference catalog. 

406 

407 Returns the candidate matched patterns and their 

408 implied rotation matrices or None. 

409 

410 Parameters 

411 ---------- 

412 src_pattern_array : `numpy.ndarray`, (N, 3) 

413 Sub selection of source 3 vectors to create a pattern from 

414 n_match : `int` 

415 Number of points to attempt to create a pattern from. Must be 

416 >= len(src_pattern_array) 

417 max_cos_theta_shift : `float` 

418 Maximum shift allowed between two patterns' centers. 

419 max_cos_rot_sq : `float` 

420 Maximum rotation between two patterns that have been shifted 

421 to have their centers on top of each other. 

422 max_dist_rad : `float` 

423 Maximum delta distance allowed between the source and reference 

424 pair distances to consider the reference pair a candidate for 

425 the source pair. Also sets the tolerance between the opening 

426 angles of the spokes when compared to the reference. 

427 

428 Return 

429 ------- 

430 output_matched_pattern : `lsst.pipe.base.Struct` 

431 Result struct with components: 

432 

433 - ``ref_candidates`` : ids of the matched pattern in the internal 

434 reference_array object (`list` of `int`). 

435 - ``src_candidates`` : Pattern ids of the sources matched 

436 (`list` of `int`). 

437 - ``shift_rot_matrix_src_to_ref`` : 3x3 matrix specifying the full 

438 shift and rotation between the reference and source objects. 

439 Rotates source into reference frame. `None` if match is not 

440 found. (`numpy.ndarray`, (3, 3)) 

441 - ``shift_rot_matrix_ref_to_src`` : 3x3 matrix specifying the full 

442 shift and rotation of the reference and source objects. Rotates 

443 reference into source frame. None if match is not found 

444 (`numpy.ndarray`, (3, 3)). 

445 - ``cos_shift`` : Magnitude of the shift found between the two 

446 patten centers. `None` if match is not found (`float`). 

447 - ``sin_rot`` : float value of the rotation to align the already 

448 shifted source pattern to the reference pattern. `None` if no match 

449 found (`float`). 

450 """ 

451 

452 # Create our place holder variables for the matched sources and 

453 # references. The source list starts with the 0th and first indexed 

454 # objects as we are guaranteed to use those and these define both 

455 # the shift and rotation of the final pattern. 

456 output_matched_pattern = pipeBase.Struct( 

457 ref_candidates=[], 

458 src_candidates=[], 

459 shift_rot_matrix=None, 

460 cos_shift=None, 

461 sin_rot=None) 

462 

463 # Create the delta vectors and distances we will need to assemble the 

464 # spokes of the pattern. 

465 src_delta_array = np.empty((len(src_pattern_array) - 1, 3)) 

466 src_delta_array[:, 0] = (src_pattern_array[1:, 0] 

467 - src_pattern_array[0, 0]) 

468 src_delta_array[:, 1] = (src_pattern_array[1:, 1] 

469 - src_pattern_array[0, 1]) 

470 src_delta_array[:, 2] = (src_pattern_array[1:, 2] 

471 - src_pattern_array[0, 2]) 

472 src_dist_array = np.sqrt(src_delta_array[:, 0]**2 

473 + src_delta_array[:, 1]**2 

474 + src_delta_array[:, 2]**2) 

475 

476 # Our first test. We search the reference dataset for pairs 

477 # that have the same length as our first source pairs to with 

478 # plus/minus the max_dist tolerance. 

479 ref_dist_index_array = self._find_candidate_reference_pairs( 

480 src_dist_array[0], self._dist_array, max_dist_rad) 

481 

482 # Start our loop over the candidate reference objects. 

483 for ref_dist_idx in ref_dist_index_array: 

484 # We have two candidates for which reference object corresponds 

485 # with the source at the center of our pattern. As such we loop 

486 # over and test both possibilities. 

487 tmp_ref_pair_list = self._id_array[ref_dist_idx] 

488 for pair_idx, ref_id in enumerate(tmp_ref_pair_list): 

489 src_candidates = [0, 1] 

490 ref_candidates = [] 

491 shift_rot_matrix = None 

492 cos_shift = None 

493 sin_rot = None 

494 # Test the angle between our candidate ref center and the 

495 # source center of our pattern. This angular distance also 

496 # defines the shift we will later use. 

497 ref_center = self._reference_array[ref_id] 

498 cos_shift = np.dot(src_pattern_array[0], ref_center) 

499 if cos_shift < max_cos_theta_shift: 

500 continue 

501 

502 # We can now append this one as a candidate. 

503 ref_candidates.append(ref_id) 

504 # Test to see which reference object to use in the pair. 

505 if pair_idx == 0: 

506 ref_candidates.append( 

507 tmp_ref_pair_list[1]) 

508 ref_delta = (self._reference_array[tmp_ref_pair_list[1]] 

509 - ref_center) 

510 else: 

511 ref_candidates.append( 

512 tmp_ref_pair_list[0]) 

513 ref_delta = (self._reference_array[tmp_ref_pair_list[0]] 

514 - ref_center) 

515 

516 # For dense fields it will be faster to compute the absolute 

517 # rotation this pair suggests first rather than saving it 

518 # after all the spokes are found. We then compute the cos^2 

519 # of the rotation and first part of the rotation matrix from 

520 # source to reference frame. 

521 test_rot_struct = self._test_rotation( 

522 src_pattern_array[0], ref_center, src_delta_array[0], 

523 ref_delta, cos_shift, max_cos_rot_sq) 

524 if test_rot_struct.cos_rot_sq is None or \ 

525 test_rot_struct.proj_ref_ctr_delta is None or \ 

526 test_rot_struct.shift_matrix is None: 

527 continue 

528 

529 # Get the data from the return struct. 

530 cos_rot_sq = test_rot_struct.cos_rot_sq 

531 proj_ref_ctr_delta = test_rot_struct.proj_ref_ctr_delta 

532 shift_matrix = test_rot_struct.shift_matrix 

533 

534 # Now that we have a candidate first spoke and reference 

535 # pattern center, we mask our future search to only those 

536 # pairs that contain our candidate reference center. 

537 tmp_ref_dist_array = self._pair_dist_array[ref_id] 

538 tmp_ref_id_array = self._pair_id_array[ref_id] 

539 

540 # Now we feed this sub data to match the spokes of 

541 # our pattern. 

542 pattern_spoke_struct = self._create_pattern_spokes( 

543 src_pattern_array[0], src_delta_array, src_dist_array, 

544 self._reference_array[ref_id], ref_id, proj_ref_ctr_delta, 

545 tmp_ref_dist_array, tmp_ref_id_array, max_dist_rad, 

546 n_match) 

547 

548 # If we don't find enough candidates we can continue to the 

549 # next reference center pair. 

550 if len(pattern_spoke_struct.ref_spoke_list) < n_match - 2 or \ 

551 len(pattern_spoke_struct.src_spoke_list) < n_match - 2: 

552 continue 

553 

554 # If we have the right number of matched ids we store these. 

555 ref_candidates.extend(pattern_spoke_struct.ref_spoke_list) 

556 src_candidates.extend(pattern_spoke_struct.src_spoke_list) 

557 

558 # We can now create our full rotation matrix for both the 

559 # shift and rotation. Reminder shift, aligns the pattern 

560 # centers, rotation rotates the spokes on top of each other. 

561 shift_rot_struct = self._create_shift_rot_matrix( 

562 cos_rot_sq, shift_matrix, src_delta_array[0], 

563 self._reference_array[ref_id], ref_delta) 

564 # If we fail to create the rotation matrix, continue to the 

565 # next objects. 

566 if shift_rot_struct.sin_rot is None or \ 

567 shift_rot_struct.shift_rot_matrix is None: 

568 continue 

569 

570 # Get the data from the return struct. 

571 sin_rot = shift_rot_struct.sin_rot 

572 shift_rot_matrix = shift_rot_struct.shift_rot_matrix 

573 

574 # Now that we have enough candidates we test to see if it 

575 # passes intermediate verify. This shifts and rotates the 

576 # source pattern into the reference frame and then tests that 

577 # each source/reference object pair is within max_dist. It also 

578 # tests the opposite rotation that is reference to source 

579 # frame. 

580 fit_shift_rot_matrix = self._intermediate_verify( 

581 src_pattern_array[src_candidates], 

582 self._reference_array[ref_candidates], 

583 shift_rot_matrix, max_dist_rad) 

584 

585 if fit_shift_rot_matrix is not None: 

586 # Fill the struct and return. 

587 output_matched_pattern.ref_candidates = ref_candidates 

588 output_matched_pattern.src_candidates = src_candidates 

589 output_matched_pattern.shift_rot_matrix = \ 

590 fit_shift_rot_matrix 

591 output_matched_pattern.cos_shift = cos_shift 

592 output_matched_pattern.sin_rot = sin_rot 

593 return output_matched_pattern 

594 

595 return output_matched_pattern 

596 

597 def _find_candidate_reference_pairs(self, src_dist, ref_dist_array, 

598 max_dist_rad): 

599 """Wrap numpy.searchsorted to find the range of reference spokes 

600 within a spoke distance tolerance of our source spoke. 

601 

602 Returns an array sorted from the smallest absolute delta distance 

603 between source and reference spoke length. This sorting increases the 

604 speed for the pattern search greatly. 

605 

606 Parameters 

607 ---------- 

608 src_dist : `float` 

609 float value of the distance we would like to search for in 

610 the reference array in radians. 

611 ref_dist_array : `numpy.ndarray`, (N,) 

612 sorted array of distances in radians. 

613 max_dist_rad : `float` 

614 maximum plus/minus search to find in the reference array in 

615 radians. 

616 

617 Return 

618 ------ 

619 tmp_diff_array : `numpy.ndarray`, (N,) 

620 indices lookup into the input ref_dist_array sorted by the 

621 difference in value to the src_dist from absolute value 

622 smallest to largest. 

623 """ 

624 # Find the index of the minimum and maximum values that satisfy 

625 # the tolerance. 

626 start_idx = np.searchsorted(ref_dist_array, src_dist - max_dist_rad) 

627 end_idx = np.searchsorted(ref_dist_array, src_dist + max_dist_rad, 

628 side='right') 

629 

630 # If these are equal there are no candidates and we exit. 

631 if start_idx == end_idx: 

632 return [] 

633 

634 # Make sure the endpoints of the input array are respected. 

635 if start_idx < 0: 

636 start_idx = 0 

637 if end_idx > ref_dist_array.shape[0]: 

638 end_idx = ref_dist_array.shape[0] 

639 

640 # Now we sort the indices from smallest absolute delta dist difference 

641 # to the largest and return the vector. This step greatly increases 

642 # the speed of the algorithm. 

643 tmp_diff_array = np.fabs(ref_dist_array[start_idx:end_idx] - src_dist) 

644 return tmp_diff_array.argsort() + start_idx 

645 

646 def _test_rotation(self, src_center, ref_center, src_delta, ref_delta, 

647 cos_shift, max_cos_rot_sq): 

648 """ Test if the rotation implied between the source 

649 pattern and reference pattern is within tolerance. To test this 

650 we need to create the first part of our spherical rotation matrix 

651 which we also return for use later. 

652 

653 Parameters 

654 ---------- 

655 src_center : `numpy.ndarray`, (N, 3) 

656 pattern. 

657 ref_center : `numpy.ndarray`, (N, 3) 

658 3 vector defining the center of the candidate reference pinwheel 

659 pattern. 

660 src_delta : `numpy.ndarray`, (N, 3) 

661 3 vector delta between the source pattern center and the end of 

662 the pinwheel spoke. 

663 ref_delta : `numpy.ndarray`, (N, 3) 

664 3 vector delta of the candidate matched reference pair 

665 cos_shift : `float` 

666 Cosine of the angle between the source and reference candidate 

667 centers. 

668 max_cos_rot_sq : `float` 

669 candidate reference pair after shifting the centers on top of each 

670 other. The function will return None if the rotation implied is 

671 greater than max_cos_rot_sq. 

672 

673 Returns 

674 ------- 

675 result : `lsst.pipe.base.Struct` 

676 Result struct with components: 

677 

678 - ``cos_rot_sq`` : magnitude of the rotation needed to align the 

679 two patterns after their centers are shifted on top of each 

680 other. `None` if rotation test fails (`float`). 

681 - ``shift_matrix`` : 3x3 rotation matrix describing the shift needed to 

682 align the source and candidate reference center. `None` if rotation 

683 test fails (`numpy.ndarray`, (N, 3)). 

684 """ 

685 

686 # Make sure the sine is a real number. 

687 if cos_shift > 1.0: 

688 cos_shift = 1. 

689 elif cos_shift < -1.0: 

690 cos_shift = -1. 

691 sin_shift = np.sqrt(1 - cos_shift ** 2) 

692 

693 # If the sine of our shift is zero we only need to use the identity 

694 # matrix for the shift. Else we construct the rotation matrix for 

695 # shift. 

696 if sin_shift > 0: 

697 rot_axis = np.cross(src_center, ref_center) 

698 rot_axis /= sin_shift 

699 shift_matrix = self._create_spherical_rotation_matrix( 

700 rot_axis, cos_shift, sin_shift) 

701 else: 

702 shift_matrix = np.identity(3) 

703 

704 # Now that we have our shift we apply it to the src delta vector 

705 # and check the rotation. 

706 rot_src_delta = np.dot(shift_matrix, src_delta) 

707 proj_src_delta = (rot_src_delta 

708 - np.dot(rot_src_delta, ref_center) * ref_center) 

709 proj_ref_delta = (ref_delta 

710 - np.dot(ref_delta, ref_center) * ref_center) 

711 cos_rot_sq = (np.dot(proj_src_delta, proj_ref_delta) ** 2 

712 / (np.dot(proj_src_delta, proj_src_delta) 

713 * np.dot(proj_ref_delta, proj_ref_delta))) 

714 # If the rotation isn't in tolerance return None. 

715 if cos_rot_sq < max_cos_rot_sq: 

716 return pipeBase.Struct( 

717 cos_rot_sq=None, 

718 proj_ref_ctr_delta=None, 

719 shift_matrix=None,) 

720 # Return the rotation angle, the plane projected reference vector, 

721 # and the first half of the full shift and rotation matrix. 

722 return pipeBase.Struct( 

723 cos_rot_sq=cos_rot_sq, 

724 proj_ref_ctr_delta=proj_ref_delta, 

725 shift_matrix=shift_matrix,) 

726 

727 def _create_spherical_rotation_matrix(self, rot_axis, cos_rotation, 

728 sin_rotion): 

729 """Construct a generalized 3D rotation matrix about a given 

730 axis. 

731 

732 Parameters 

733 ---------- 

734 rot_axis : `numpy.ndarray`, (3,) 

735 3 vector defining the axis to rotate about. 

736 cos_rotation : `float` 

737 cosine of the rotation angle. 

738 sin_rotation : `float` 

739 sine of the rotation angle. 

740 

741 Return 

742 ------ 

743 shift_matrix : `numpy.ndarray`, (3, 3) 

744 3x3 spherical, rotation matrix. 

745 """ 

746 

747 rot_cross_matrix = np.array( 

748 [[0., -rot_axis[2], rot_axis[1]], 

749 [rot_axis[2], 0., -rot_axis[0]], 

750 [-rot_axis[1], rot_axis[0], 0.]], dtype=np.float64) 

751 shift_matrix = (cos_rotation*np.identity(3) 

752 + sin_rotion*rot_cross_matrix 

753 + (1. - cos_rotation)*np.outer(rot_axis, rot_axis)) 

754 

755 return shift_matrix 

756 

757 def _create_pattern_spokes(self, src_ctr, src_delta_array, src_dist_array, 

758 ref_ctr, ref_ctr_id, proj_ref_ctr_delta, 

759 ref_dist_array, ref_id_array, max_dist_rad, 

760 n_match): 

761 """ Create the individual spokes that make up the pattern now that the 

762 shift and rotation are within tolerance. 

763 

764 If we can't create a valid pattern we exit early. 

765 

766 Parameters 

767 ---------- 

768 src_ctr : `numpy.ndarray`, (3,) 

769 3 vector of the source pinwheel center 

770 src_delta_array : `numpy.ndarray`, (N, 3) 

771 Array of 3 vector deltas between the source center and the pairs 

772 that make up the remaining spokes of the pinwheel 

773 src_dist_array : `numpy.ndarray`, (N, 3) 

774 Array of the distances of each src_delta in the pinwheel 

775 ref_ctr : `numpy.ndarray`, (3,) 

776 3 vector of the candidate reference center 

777 ref_ctr_id : `int` 

778 id of the ref_ctr in the master reference array 

779 proj_ref_ctr_delta : `numpy.ndarray`, (3,) 

780 Plane projected 3 vector formed from the center point of the 

781 candidate pin-wheel and the second point in the pattern to create 

782 the first spoke pair. This is the candidate pair that was matched 

783 in the main _construct_pattern_and_shift_rot_matrix loop 

784 ref_dist_array : `numpy.ndarray`, (N,) 

785 Array of vector distances for each of the reference pairs 

786 ref_id_array : `numpy.ndarray`, (N,) 

787 Array of id lookups into the master reference array that our 

788 center id object is paired with. 

789 max_dist_rad : `float` 

790 Maximum search distance 

791 n_match : `int` 

792 Number of source deltas that must be matched into the reference 

793 deltas in order to consider this a successful pattern match. 

794 

795 Returns 

796 ------- 

797 output_spokes : `lsst.pipe.base.Struct` 

798 Result struct with components: 

799 

800 - ``ref_spoke_list`` : list of ints specifying ids into the master 

801 reference array (`list` of `int`). 

802 - ``src_spoke_list`` : list of ints specifying indices into the 

803 current source pattern that is being tested (`list` of `int`). 

804 """ 

805 # Struct where we will be putting our results. 

806 output_spokes = pipeBase.Struct( 

807 ref_spoke_list=[], 

808 src_spoke_list=[],) 

809 

810 # Counter for number of spokes we failed to find a reference 

811 # candidate for. We break the loop if we haven't found enough. 

812 n_fail = 0 

813 ref_spoke_list = [] 

814 src_spoke_list = [] 

815 

816 # Plane project the center/first spoke of the source pattern using 

817 # the center vector of the pattern as normal. 

818 proj_src_ctr_delta = (src_delta_array[0] 

819 - np.dot(src_delta_array[0], src_ctr) * src_ctr) 

820 proj_src_ctr_dist_sq = np.dot(proj_src_ctr_delta, proj_src_ctr_delta) 

821 

822 # Pre-compute the squared length of the projected reference vector. 

823 proj_ref_ctr_dist_sq = np.dot(proj_ref_ctr_delta, proj_ref_ctr_delta) 

824 

825 # Loop over the source pairs. 

826 for src_idx in range(1, len(src_dist_array)): 

827 if n_fail > len(src_dist_array) - (n_match - 1): 

828 break 

829 

830 # Given our length tolerance we can use it to compute a tolerance 

831 # on the angle between our spoke. 

832 src_sin_tol = (max_dist_rad 

833 / (src_dist_array[src_idx] + max_dist_rad)) 

834 

835 # Test if the small angle approximation will still hold. This is 

836 # defined as when sin(theta) ~= theta to within 0.1% of each 

837 # other. If the implied opening angle is too large we set it to 

838 # the 0.1% threshold. 

839 max_sin_tol = 0.0447 

840 if src_sin_tol > max_sin_tol: 

841 src_sin_tol = max_sin_tol 

842 

843 # Plane project the candidate source spoke and compute the cosine 

844 # and sine of the opening angle. 

845 proj_src_delta = ( 

846 src_delta_array[src_idx] 

847 - np.dot(src_delta_array[src_idx], src_ctr) * src_ctr) 

848 geom_dist_src = np.sqrt( 

849 np.dot(proj_src_delta, proj_src_delta) 

850 * proj_src_ctr_dist_sq) 

851 

852 # Compute cosine and sine of the delta vector opening angle. 

853 cos_theta_src = (np.dot(proj_src_delta, proj_src_ctr_delta) 

854 / geom_dist_src) 

855 cross_src = (np.cross(proj_src_delta, proj_src_ctr_delta) 

856 / geom_dist_src) 

857 sin_theta_src = np.dot(cross_src, src_ctr) 

858 

859 # Find the reference pairs that include our candidate pattern 

860 # center and sort them in increasing delta 

861 ref_dist_idx_array = self._find_candidate_reference_pairs( 

862 src_dist_array[src_idx], ref_dist_array, max_dist_rad) 

863 

864 # Test the spokes and return the id of the reference object. 

865 # Return None if no match is found. 

866 ref_id = self._test_spoke( 

867 cos_theta_src, 

868 sin_theta_src, 

869 ref_ctr, 

870 ref_ctr_id, 

871 proj_ref_ctr_delta, 

872 proj_ref_ctr_dist_sq, 

873 ref_dist_idx_array, 

874 ref_id_array, 

875 src_sin_tol) 

876 if ref_id is None: 

877 n_fail += 1 

878 continue 

879 

880 # Append the successful indices to our list. The src_idx needs 

881 # an extra iteration to skip the first and second source objects. 

882 ref_spoke_list.append(ref_id) 

883 src_spoke_list.append(src_idx + 1) 

884 # If we found enough reference objects we can return early. This is 

885 # n_match - 2 as we already have 2 source objects matched into the 

886 # reference data. 

887 if len(ref_spoke_list) >= n_match - 2: 

888 # Set the struct data and return the struct. 

889 output_spokes.ref_spoke_list = ref_spoke_list 

890 output_spokes.src_spoke_list = src_spoke_list 

891 return output_spokes 

892 

893 return output_spokes 

894 

895 def _test_spoke(self, cos_theta_src, sin_theta_src, ref_ctr, ref_ctr_id, 

896 proj_ref_ctr_delta, proj_ref_ctr_dist_sq, 

897 ref_dist_idx_array, ref_id_array, src_sin_tol): 

898 """Test the opening angle between the first spoke of our pattern 

899 for the source object against the reference object. 

900 

901 This method makes heavy use of the small angle approximation to perform 

902 the comparison. 

903 

904 Parameters 

905 ---------- 

906 cos_theta_src : `float` 

907 Cosine of the angle between the current candidate source spoke and 

908 the first spoke. 

909 sin_theta_src : `float` 

910 Sine of the angle between the current candidate source spoke and 

911 the first spoke. 

912 ref_ctr : `numpy.ndarray`, (3,) 

913 3 vector of the candidate reference center 

914 ref_ctr_id : `int` 

915 id lookup of the ref_ctr into the master reference array 

916 proj_ref_ctr_delta : `float` 

917 Plane projected first spoke in the reference pattern using the 

918 pattern center as normal. 

919 proj_ref_ctr_dist_sq : `float` 

920 Squared length of the projected vector. 

921 ref_dist_idx_array : `numpy.ndarray`, (N,) 

922 Indices sorted by the delta distance between the source 

923 spoke we are trying to test and the candidate reference 

924 spokes. 

925 ref_id_array : `numpy.ndarray`, (N,) 

926 Array of id lookups into the master reference array that our 

927 center id object is paired with. 

928 src_sin_tol : `float` 

929 Sine of tolerance allowed between source and reference spoke 

930 opening angles. 

931 

932 Returns 

933 ------- 

934 id : `int` 

935 If we can not find a candidate spoke we return `None` else we 

936 return an int id into the master reference array. 

937 """ 

938 

939 # Loop over our candidate reference objects. 

940 for ref_dist_idx in ref_dist_idx_array: 

941 # Compute the delta vector from the pattern center. 

942 ref_delta = (self._reference_array[ref_id_array[ref_dist_idx]] 

943 - ref_ctr) 

944 # Compute the cos between our "center" reference vector and the 

945 # current reference candidate. 

946 proj_ref_delta = ref_delta - np.dot(ref_delta, ref_ctr) * ref_ctr 

947 geom_dist_ref = np.sqrt(proj_ref_ctr_dist_sq 

948 * np.dot(proj_ref_delta, proj_ref_delta)) 

949 cos_theta_ref = (np.dot(proj_ref_delta, proj_ref_ctr_delta) 

950 / geom_dist_ref) 

951 

952 # Make sure we can safely make the comparison in case 

953 # our "center" and candidate vectors are mostly aligned. 

954 if cos_theta_ref ** 2 < (1 - src_sin_tol ** 2): 

955 cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2 

956 / (1 - cos_theta_ref ** 2)) 

957 else: 

958 cos_sq_comparison = ((cos_theta_src - cos_theta_ref) ** 2 

959 / src_sin_tol ** 2) 

960 # Test the difference of the cosine of the reference angle against 

961 # the source angle. Assumes that the delta between the two is 

962 # small. 

963 if cos_sq_comparison > src_sin_tol ** 2: 

964 continue 

965 

966 # The cosine tests the magnitude of the angle but not 

967 # its direction. To do that we need to know the sine as well. 

968 # This cross product calculation does that. 

969 cross_ref = (np.cross(proj_ref_delta, proj_ref_ctr_delta) 

970 / geom_dist_ref) 

971 sin_theta_ref = np.dot(cross_ref, ref_ctr) 

972 

973 # Check the value of the cos again to make sure that it is not 

974 # near zero. 

975 if abs(cos_theta_src) < src_sin_tol: 

976 sin_comparison = (sin_theta_src - sin_theta_ref) / src_sin_tol 

977 else: 

978 sin_comparison = \ 

979 (sin_theta_src - sin_theta_ref) / cos_theta_ref 

980 

981 if abs(sin_comparison) > src_sin_tol: 

982 continue 

983 

984 # Return the correct id of the candidate we found. 

985 return ref_id_array[ref_dist_idx] 

986 

987 return None 

988 

989 def _create_shift_rot_matrix(self, cos_rot_sq, shift_matrix, src_delta, 

990 ref_ctr, ref_delta): 

991 """ Create the final part of our spherical rotation matrix. 

992 

993 Parameters 

994 ---------- 

995 cos_rot_sq : `float` 

996 cosine of the rotation needed to align our source and reference 

997 candidate patterns. 

998 shift_matrix : `numpy.ndarray`, (3, 3) 

999 3x3 rotation matrix for shifting the source pattern center on top 

1000 of the candidate reference pattern center. 

1001 src_delta : `numpy.ndarray`, (3,) 

1002 3 vector delta of representing the first spoke of the source 

1003 pattern 

1004 ref_ctr : `numpy.ndarray`, (3,) 

1005 3 vector on the unit-sphere representing the center of our 

1006 reference pattern. 

1007 ref_delta : `numpy.ndarray`, (3,) 

1008 3 vector delta made by the first pair of the reference pattern. 

1009 

1010 Returns 

1011 ------- 

1012 result : `lsst.pipe.base.Struct` 

1013 Result struct with components: 

1014 

1015 - ``sin_rot`` : float sine of the amount of rotation between the 

1016 source and reference pattern. We use sine here as it is 

1017 signed and tells us the chirality of the rotation (`float`). 

1018 - ``shift_rot_matrix`` : float array representing the 3x3 rotation 

1019 matrix that takes the source pattern and shifts and rotates 

1020 it to align with the reference pattern (`numpy.ndarray`, (3,3)). 

1021 """ 

1022 cos_rot = np.sqrt(cos_rot_sq) 

1023 rot_src_delta = np.dot(shift_matrix, src_delta) 

1024 delta_dot_cross = np.dot(np.cross(rot_src_delta, ref_delta), ref_ctr) 

1025 

1026 sin_rot = np.sign(delta_dot_cross) * np.sqrt(1 - cos_rot_sq) 

1027 rot_matrix = self._create_spherical_rotation_matrix( 

1028 ref_ctr, cos_rot, sin_rot) 

1029 

1030 shift_rot_matrix = np.dot(rot_matrix, shift_matrix) 

1031 

1032 return pipeBase.Struct( 

1033 sin_rot=sin_rot, 

1034 shift_rot_matrix=shift_rot_matrix,) 

1035 

1036 def _intermediate_verify(self, src_pattern, ref_pattern, shift_rot_matrix, 

1037 max_dist_rad): 

1038 """ Perform an intermediate verify step. 

1039 

1040 Rotate the matches references into the source frame and test their 

1041 distances against tolerance. Only return true if all points are within 

1042 tolerance. 

1043 

1044 Parameters 

1045 ---------- 

1046 src_pattern : `numpy.ndarray`, (N,3) 

1047 Array of 3 vectors representing the points that make up our source 

1048 pinwheel pattern. 

1049 ref_pattern : `numpy.ndarray`, (N,3) 

1050 Array of 3 vectors representing our candidate reference pinwheel 

1051 pattern. 

1052 shift_rot_matrix : `numpy.ndarray`, (3,3) 

1053 3x3 rotation matrix that takes the source objects and rotates them 

1054 onto the frame of the reference objects 

1055 max_dist_rad : `float` 

1056 Maximum distance allowed to consider two objects the same. 

1057 

1058 Returns 

1059 ------- 

1060 fit_shift_rot_matrix : `numpy.ndarray`, (3,3) 

1061 Return the fitted shift/rotation matrix if all of the points in our 

1062 source pattern are within max_dist_rad of their matched reference 

1063 objects. Returns None if this criteria is not satisfied. 

1064 """ 

1065 if len(src_pattern) != len(ref_pattern): 

1066 raise ValueError( 

1067 "Source pattern length does not match ref pattern.\n" 

1068 "\t source pattern len=%i, reference pattern len=%i" % 

1069 (len(src_pattern), len(ref_pattern))) 

1070 

1071 if self._intermediate_verify_comparison( 

1072 src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad): 

1073 # Now that we know our initial shift and rot matrix is valid we 

1074 # want to fit the implied transform using all points from 

1075 # our pattern. This is a more robust rotation matrix as our 

1076 # initial matrix only used the first 2 points from the source 

1077 # pattern to estimate the shift and rotation. The matrix we fit 

1078 # are allowed to be non unitary but need to preserve the length of 

1079 # the rotated vectors to within the match tolerance. 

1080 fit_shift_rot_matrix = least_squares( 

1081 _rotation_matrix_chi_sq, 

1082 x0=shift_rot_matrix.flatten(), 

1083 args=(src_pattern, ref_pattern, max_dist_rad) 

1084 ).x.reshape((3, 3)) 

1085 # Do another verify in case the fits have wondered off. 

1086 if self._intermediate_verify_comparison( 

1087 src_pattern, ref_pattern, fit_shift_rot_matrix, 

1088 max_dist_rad): 

1089 return fit_shift_rot_matrix 

1090 

1091 return None 

1092 

1093 def _intermediate_verify_comparison(self, pattern_a, pattern_b, 

1094 shift_rot_matrix, max_dist_rad): 

1095 """Test the input rotation matrix against one input pattern and 

1096 a second one. 

1097 

1098 If every point in the pattern after rotation is within a distance of 

1099 max_dist_rad to its candidate point in the other pattern, we return 

1100 True. 

1101 

1102 Parameters 

1103 ---------- 

1104 pattern_a : `numpy.ndarray`, (N,3) 

1105 Array of 3 vectors representing the points that make up our source 

1106 pinwheel pattern. 

1107 pattern_b : `numpy.ndarray`, (N,3) 

1108 Array of 3 vectors representing our candidate reference pinwheel 

1109 pattern. 

1110 shift_rot_matrix : `numpy.ndarray`, (3,3) 

1111 3x3 rotation matrix that takes the source objects and rotates them 

1112 onto the frame of the reference objects 

1113 max_dist_rad : `float` 

1114 Maximum distance allowed to consider two objects the same. 

1115 

1116 

1117 Returns 

1118 ------- 

1119 bool 

1120 True if all rotated source points are within max_dist_rad of 

1121 the candidate references matches. 

1122 """ 

1123 shifted_pattern_a = np.dot(shift_rot_matrix, 

1124 pattern_a.transpose()).transpose() 

1125 tmp_delta_array = shifted_pattern_a - pattern_b 

1126 tmp_dist_array = (tmp_delta_array[:, 0] ** 2 

1127 + tmp_delta_array[:, 1] ** 2 

1128 + tmp_delta_array[:, 2] ** 2) 

1129 return np.all(tmp_dist_array < max_dist_rad ** 2) 

1130 

1131 def _test_pattern_lengths(self, test_pattern, max_dist_rad): 

1132 """ Test that the all vectors in a pattern are unit length within 

1133 tolerance. 

1134 

1135 This is useful for assuring the non unitary transforms do not contain 

1136 too much distortion. 

1137 

1138 Parameters 

1139 ---------- 

1140 test_pattern : `numpy.ndarray`, (N, 3) 

1141 Test vectors at the maximum and minimum x, y, z extents. 

1142 max_dist_rad : `float` 

1143 maximum distance in radians to consider two points "agreeing" on 

1144 a rotation 

1145 

1146 Returns 

1147 ------- 

1148 test : `bool` 

1149 Length tests pass. 

1150 """ 

1151 dists = (test_pattern[:, 0] ** 2 

1152 + test_pattern[:, 1] ** 2 

1153 + test_pattern[:, 2] ** 2) 

1154 return np.all( 

1155 np.logical_and((1 - max_dist_rad) ** 2 < dists, 

1156 dists < (1 + max_dist_rad) ** 2)) 

1157 

1158 def _test_rotation_agreement(self, rot_vects, max_dist_rad): 

1159 """ Test this rotation against the previous N found and return 

1160 the number that a agree within tolerance to where our test 

1161 points are. 

1162 

1163 Parameters 

1164 ---------- 

1165 rot_vects : `numpy.ndarray`, (N, 3) 

1166 Arrays of rotated 3 vectors representing the maximum x, y, 

1167 z extent on the unit sphere of the input source objects rotated by 

1168 the candidate rotations into the reference frame. 

1169 max_dist_rad : `float` 

1170 maximum distance in radians to consider two points "agreeing" on 

1171 a rotation 

1172 

1173 Returns 

1174 ------- 

1175 tot_consent : `int` 

1176 Number of candidate rotations that agree for all of the rotated 

1177 test 3 vectors. 

1178 """ 

1179 

1180 self.log.debug("Comparing pattern %i to previous %i rotations..." % 

1181 (rot_vects[-1][-1], len(rot_vects) - 1)) 

1182 

1183 tot_consent = 0 

1184 for rot_idx in range(max((len(rot_vects) - 1), 0)): 

1185 tmp_dist_list = [] 

1186 for vect_idx in range(len(rot_vects[rot_idx]) - 1): 

1187 tmp_delta_vect = (rot_vects[rot_idx][vect_idx] 

1188 - rot_vects[-1][vect_idx]) 

1189 tmp_dist_list.append( 

1190 np.dot(tmp_delta_vect, tmp_delta_vect)) 

1191 if np.all(np.array(tmp_dist_list) < max_dist_rad ** 2): 

1192 tot_consent += 1 

1193 return tot_consent 

1194 

1195 def _final_verify(self, 

1196 source_array, 

1197 shift_rot_matrix, 

1198 max_dist_rad, 

1199 min_matches): 

1200 """Match the all sources into the reference catalog using the shift/rot 

1201 matrix. 

1202 

1203 After the initial shift/rot matrix is found, we refit the shift/rot 

1204 matrix using the matches the initial matrix produces to find a more 

1205 stable solution. 

1206 

1207 Parameters 

1208 ---------- 

1209 source_array : `numpy.ndarray` (N, 3) 

1210 3-vector positions on the unit-sphere representing the sources to 

1211 match 

1212 shift_rot_matrix : `numpy.ndarray` (3, 3) 

1213 Rotation matrix representing inferred shift/rotation of the 

1214 sources onto the reference catalog. Matrix need not be unitary. 

1215 max_dist_rad : `float` 

1216 Maximum distance allowed for a match. 

1217 min_matches : `int` 

1218 Minimum number of matched objects required to consider the 

1219 match good. 

1220 

1221 Returns 

1222 ------- 

1223 output_struct : `lsst.pipe.base.Struct` 

1224 Result struct with components: 

1225 

1226 - ``match_ids`` : Pairs of indexes into the source and reference 

1227 data respectively defining a match (`numpy.ndarray`, (N, 2)). 

1228 - ``distances_rad`` : distances to between the matched objects in 

1229 the shift/rotated frame. (`numpy.ndarray`, (N,)). 

1230 - ``max_dist_rad`` : Value of the max matched distance. Either 

1231 returning the input value of the 2 sigma clipped value of the 

1232 shift/rotated distances. (`float`). 

1233 """ 

1234 output_struct = pipeBase.Struct( 

1235 match_ids=None, 

1236 distances_rad=None, 

1237 max_dist_rad=None, 

1238 ) 

1239 

1240 # Perform an iterative final verify. 

1241 match_sources_struct = self._match_sources(source_array, 

1242 shift_rot_matrix) 

1243 cut_ids = match_sources_struct.match_ids[ 

1244 match_sources_struct.distances_rad < max_dist_rad] 

1245 

1246 n_matched = len(cut_ids) 

1247 clipped_struct = self._clip_distances( 

1248 match_sources_struct.distances_rad) 

1249 n_matched_clipped = clipped_struct.n_matched_clipped 

1250 

1251 if n_matched < min_matches or n_matched_clipped < min_matches: 

1252 return output_struct 

1253 

1254 # The shift/rotation matrix returned by 

1255 # ``_construct_pattern_and_shift_rot_matrix``, above, was 

1256 # based on only six points. Here, we refine that result by 

1257 # using all of the good matches from the “final verification” 

1258 # step, above. This will produce a more consistent result. 

1259 fit_shift_rot_matrix = least_squares( 

1260 _rotation_matrix_chi_sq, 

1261 x0=shift_rot_matrix.flatten(), 

1262 args=(source_array[cut_ids[:, 0], :3], 

1263 self._reference_array[cut_ids[:, 1], :3], 

1264 max_dist_rad) 

1265 ).x.reshape((3, 3)) 

1266 

1267 # Redo the matching using the newly fit shift/rotation matrix. 

1268 match_sources_struct = self._match_sources( 

1269 source_array, fit_shift_rot_matrix) 

1270 

1271 # Double check the match distances to make sure enough matches 

1272 # survive still. We'll just overwrite the previously used variables. 

1273 n_matched = np.sum( 

1274 match_sources_struct.distances_rad < max_dist_rad) 

1275 clipped_struct = self._clip_distances( 

1276 match_sources_struct.distances_rad) 

1277 n_matched_clipped = clipped_struct.n_matched_clipped 

1278 clipped_max_dist = clipped_struct.clipped_max_dist 

1279 

1280 if n_matched < min_matches or n_matched_clipped < min_matches: 

1281 return output_struct 

1282 

1283 # Store our matches in the output struct. Decide between the clipped 

1284 # distance and the input max dist based on which is smaller. 

1285 output_struct.match_ids = match_sources_struct.match_ids 

1286 output_struct.distances_rad = match_sources_struct.distances_rad 

1287 if clipped_max_dist < max_dist_rad: 

1288 output_struct.max_dist_rad = clipped_max_dist 

1289 else: 

1290 output_struct.max_dist_rad = max_dist_rad 

1291 

1292 return output_struct 

1293 

1294 def _clip_distances(self, distances_rad): 

1295 """Compute a clipped max distance and calculate the number of pairs 

1296 that pass the clipped dist. 

1297 

1298 Parameters 

1299 ---------- 

1300 distances_rad : `numpy.ndarray`, (N,) 

1301 Distances between pairs. 

1302 

1303 Returns 

1304 ------- 

1305 output_struct : `lsst.pipe.base.Struct` 

1306 Result struct with components: 

1307 

1308 - ``n_matched_clipped`` : Number of pairs that survive the 

1309 clipping on distance. (`float`) 

1310 - ``clipped_max_dist`` : Maximum distance after clipping. 

1311 (`float`). 

1312 """ 

1313 clipped_dists, _, clipped_max_dist = sigmaclip( 

1314 distances_rad, 

1315 low=100, 

1316 high=2) 

1317 # Check clipped distances. The minimum value here 

1318 # prevents over convergence on perfect test data. 

1319 if clipped_max_dist < 1e-16: 

1320 clipped_max_dist = 1e-16 

1321 n_matched_clipped = np.sum(distances_rad < clipped_max_dist) 

1322 else: 

1323 n_matched_clipped = len(clipped_dists) 

1324 

1325 return pipeBase.Struct(n_matched_clipped=n_matched_clipped, 

1326 clipped_max_dist=clipped_max_dist) 

1327 

1328 def _match_sources(self, 

1329 source_array, 

1330 shift_rot_matrix): 

1331 """ Shift both the reference and source catalog to the the respective 

1332 frames and find their nearest neighbor using a kdTree. 

1333 

1334 Removes all matches who do not agree when either the reference or 

1335 source catalog is rotated. Cuts on a maximum distance are left to an 

1336 external function. 

1337 

1338 Parameters 

1339 ---------- 

1340 source_array : `numpy.ndarray`, (N, 3) 

1341 array of 3 vectors representing the source objects we are trying 

1342 to match into the source catalog. 

1343 shift_rot_matrix : `numpy.ndarray`, (3, 3) 

1344 3x3 rotation matrix that performs the spherical rotation from the 

1345 source frame into the reference frame. 

1346 

1347 Returns 

1348 ------- 

1349 results : `lsst.pipe.base.Struct` 

1350 Result struct with components: 

1351 

1352 - ``matches`` : array of integer ids into the source and 

1353 reference arrays. Matches are only returned for those that 

1354 satisfy the distance and handshake criteria 

1355 (`numpy.ndarray`, (N, 2)). 

1356 - ``distances`` : Distances between each match in radians after 

1357 the shift and rotation is applied (`numpy.ndarray`, (N)). 

1358 """ 

1359 shifted_references = np.dot( 

1360 np.linalg.inv(shift_rot_matrix), 

1361 self._reference_array.transpose()).transpose() 

1362 shifted_sources = np.dot( 

1363 shift_rot_matrix, 

1364 source_array.transpose()).transpose() 

1365 

1366 ref_matches = np.empty((len(shifted_references), 2), 

1367 dtype=np.uint16) 

1368 src_matches = np.empty((len(shifted_sources), 2), 

1369 dtype=np.uint16) 

1370 

1371 ref_matches[:, 1] = np.arange(len(shifted_references), 

1372 dtype=np.uint16) 

1373 src_matches[:, 0] = np.arange(len(shifted_sources), 

1374 dtype=np.uint16) 

1375 

1376 ref_kdtree = cKDTree(self._reference_array) 

1377 src_kdtree = cKDTree(source_array) 

1378 

1379 ref_to_src_dist, tmp_ref_to_src_idx = \ 

1380 src_kdtree.query(shifted_references) 

1381 src_to_ref_dist, tmp_src_to_ref_idx = \ 

1382 ref_kdtree.query(shifted_sources) 

1383 

1384 ref_matches[:, 0] = tmp_ref_to_src_idx 

1385 src_matches[:, 1] = tmp_src_to_ref_idx 

1386 

1387 handshake_mask = self._handshake_match(src_matches, ref_matches) 

1388 return pipeBase.Struct( 

1389 match_ids=src_matches[handshake_mask], 

1390 distances_rad=src_to_ref_dist[handshake_mask],) 

1391 

1392 def _handshake_match(self, matches_src, matches_ref): 

1393 """Return only those matches where both the source 

1394 and reference objects agree they they are each others' 

1395 nearest neighbor. 

1396 

1397 Parameters 

1398 ---------- 

1399 matches_src : `numpy.ndarray`, (N, 2) 

1400 int array of nearest neighbor matches between shifted and 

1401 rotated reference objects matched into the sources. 

1402 matches_ref : `numpy.ndarray`, (N, 2) 

1403 int array of nearest neighbor matches between shifted and 

1404 rotated source objects matched into the references. 

1405 Return 

1406 ------ 

1407 handshake_mask_array : `numpy.ndarray`, (N,) 

1408 Return the array positions where the two match catalogs agree. 

1409 """ 

1410 handshake_mask_array = np.zeros(len(matches_src), dtype=np.bool) 

1411 

1412 for src_match_idx, match in enumerate(matches_src): 

1413 ref_match_idx = np.searchsorted(matches_ref[:, 1], match[1]) 

1414 if match[0] == matches_ref[ref_match_idx, 0]: 

1415 handshake_mask_array[src_match_idx] = True 

1416 return handshake_mask_array