lsst.meas.astrom  22.0.1-11-g21ef37a+cb53701f20
pessimisticPatternMatcherUtils.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 
3 /*
4  * This file is part of meas_astrom.
5  *
6  * Developed for the LSST Data Management System.
7  * This product includes software developed by the LSST Project
8  * (https://www.lsst.org).
9  * See the COPYRIGHT file at the top-level directory of this distribution
10  * for details of code ownership.
11  *
12  * This program is free software: you can redistribute it and/or modify
13  * it under the terms of the GNU General Public License as published by
14  * the Free Software Foundation, either version 3 of the License, or
15  * (at your option) any later version.
16  *
17  * This program is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  * GNU General Public License for more details.
21  *
22  * You should have received a copy of the GNU General Public License
23  * along with this program. If not, see <https://www.gnu.org/licenses/>.
24  */
25 
26 #include <cmath>
27 #include <Eigen/Dense>
28 #include "ndarray/eigen.h"
30 
31 namespace lsst {
32 namespace meas {
33 namespace astrom {
34 
35 int check_spoke(double cos_theta_src, double sin_theta_src, ndarray::Array<double, 1, 1> const& ref_ctr,
36  ndarray::Array<double, 1, 1> const& proj_ref_ctr_delta, double proj_ref_ctr_dist_sq,
37  ndarray::Array<long int, 1, 1> const& ref_dist_idx_array,
38  ndarray::Array<uint16_t, 1, 1> const& ref_id_array,
39  ndarray::Array<double, 2, 1> const& reference_array, double src_sin_tol) {
40  // Loop over our candidate reference objects. ref_dist_idx_array is a view into ref_id_array
41  // and is not the same length as ref_id_array.
42  for (auto idx = ref_dist_idx_array.begin(); idx != ref_dist_idx_array.end(); idx++) {
43  // Compute the delta vector from the pattern center.
44  unsigned int ref_id = ref_id_array[*idx];
45  ndarray::Array<double, 1, 1> ref_delta = copy(reference_array[ref_id] - ref_ctr);
46 
47  double ref_dot = ndarray::asEigenMatrix(ref_delta).dot(ndarray::asEigenMatrix(ref_ctr));
48  ndarray::Array<double, 1, 1> proj_ref_delta = copy(ref_delta - ref_dot * ref_ctr);
49  // Compute the cos between our "center" reference vector and the
50  // current reference candidate.
51  double proj_delta_dist_sq =
52  ndarray::asEigenMatrix(proj_ref_delta).dot(ndarray::asEigenMatrix(proj_ref_delta));
53  double geom_dist_ref = sqrt(proj_ref_ctr_dist_sq * proj_delta_dist_sq);
54  double cos_theta_ref =
55  ndarray::asEigenMatrix(proj_ref_delta).dot(ndarray::asEigenMatrix(proj_ref_ctr_delta)) /
56  geom_dist_ref;
57 
58  // Make sure we can safely make the comparison in case
59  // our "center" and candidate vectors are mostly aligned.
60  double cos_sq_comparison;
61  if (cos_theta_ref * cos_theta_ref < (1 - src_sin_tol * src_sin_tol)) {
62  cos_sq_comparison = (cos_theta_src - cos_theta_ref) * (cos_theta_src - cos_theta_ref) /
63  (1 - cos_theta_ref * cos_theta_ref);
64  } else {
65  cos_sq_comparison = (cos_theta_src - cos_theta_ref) * (cos_theta_src - cos_theta_ref) /
66  (src_sin_tol * src_sin_tol);
67  }
68  // Test the difference of the cosine of the reference angle against
69  // the source angle. Assumes that the delta between the two is
70  // small.
71  if (cos_sq_comparison > src_sin_tol * src_sin_tol) {
72  continue;
73  }
74  // The cosine tests the magnitude of the angle but not
75  // its direction. To do that we need to know the sine as well.
76  // This cross product calculation does that.
77  Eigen::Vector3d cross_ref = ndarray::asEigenMatrix(proj_ref_delta)
78  .head<3>()
79  .cross(ndarray::asEigenMatrix(proj_ref_ctr_delta).head<3>()) /
80  geom_dist_ref;
81  double sin_theta_ref = cross_ref.dot(ndarray::asEigenMatrix(ref_ctr));
82  // Check the value of the cos again to make sure that it is not
83  // near zero.
84  double sin_comparison;
85  if (abs(cos_theta_src) < src_sin_tol) {
86  sin_comparison = (sin_theta_src - sin_theta_ref) / src_sin_tol;
87  } else {
88  sin_comparison = (sin_theta_src - sin_theta_ref) / cos_theta_ref;
89  }
90 
91  // Return the correct id of the candidate we found.
92  if (abs(sin_comparison) < src_sin_tol) {
93  return ref_id;
94  }
95  }
96  return -1;
97 }
98 
99 } // namespace astrom
100 } // namespace meas
101 } // namespace lsst
int check_spoke(double cos_theta_src, double sin_theta_src, ndarray::Array< double, 1, 1 > const &ref_ctr, ndarray::Array< double, 1, 1 > const &proj_ref_ctr_delta, double proj_ref_ctr_dist_sq, ndarray::Array< long int, 1, 1 > const &ref_dist_idx_array, ndarray::Array< uint16_t, 1, 1 > const &ref_id_array, ndarray::Array< double, 2, 1 > const &reference_array, double src_sin_tol)
Check the opening angle between the first spoke of our pattern for the source object against the refe...