LSST Applications g0f08755f38+9c285cab97,g1635faa6d4+13f3999e92,g1653933729+a8ce1bb630,g1a0ca8cf93+bf6eb00ceb,g28da252d5a+0829b12dee,g29321ee8c0+5700dc9eac,g2bbee38e9b+9634bc57db,g2bc492864f+9634bc57db,g2cdde0e794+c2c89b37c4,g3156d2b45e+41e33cbcdc,g347aa1857d+9634bc57db,g35bb328faa+a8ce1bb630,g3a166c0a6a+9634bc57db,g3e281a1b8c+9f2c4e2fc3,g414038480c+077ccc18e7,g41af890bb2+fde0dd39b6,g5fbc88fb19+17cd334064,g781aacb6e4+a8ce1bb630,g80478fca09+55a9465950,g82479be7b0+d730eedb7d,g858d7b2824+9c285cab97,g9125e01d80+a8ce1bb630,g9726552aa6+10f999ec6a,ga5288a1d22+2a84bb7594,gacf8899fa4+c69c5206e8,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gc28159a63d+9634bc57db,gcf0d15dbbd+4b7d09cae4,gda3e153d99+9c285cab97,gda6a2b7d83+4b7d09cae4,gdaeeff99f8+1711a396fd,ge2409df99d+5e831397f4,ge79ae78c31+9634bc57db,gf0baf85859+147a0692ba,gf3967379c6+41c94011de,gf3fb38a9a8+8f07a9901b,gfb92a5be7c+9c285cab97,w.2024.46
LSST Data Management Base Package
Loading...
Searching...
No Matches
Namespaces | Classes | Typedefs | Functions
lsst::meas::astrom Namespace Reference

Namespaces

namespace  approximateWcs
 
namespace  astrometry
 
namespace  denormalizeMatches
 
namespace  detail
 
namespace  directMatch
 
namespace  display
 
namespace  exceptions
 
namespace  fitAffineWcs
 
namespace  fitSipDistortion
 
namespace  fitTanSipWcs
 
namespace  match_probabilistic_task
 
namespace  matcher_probabilistic
 
namespace  matchOptimisticB
 
namespace  matchOptimisticBTask
 
namespace  matchPessimisticB
 
namespace  pessimistic_pattern_matcher_b_3D
 
namespace  ref_match
 
namespace  setMatchDistance
 
namespace  sip
 
namespace  verifyWcs
 
namespace  version
 

Classes

struct  MatchOptimisticBControl
 
class  OutlierRejectionControl
 Control object for outlier rejection in ScaledPolynomialTransformFitter. More...
 
struct  PatternResult
 Result of construct_pattern_and_shift_rot_matrix(), containing the matched sources, rotation matrix, and a success flag. More...
 
class  PolynomialTransform
 A 2-d coordinate transform represented by a pair of standard polynomials (one for each coordinate). More...
 
struct  ProxyPair
 
struct  RecordProxy
 A wrapper around a SimpleRecord or SourceRecord that allows us to record a pixel position in a way that is independent of the record type. More...
 
struct  RotationTestResult
 Result of test_rotation(), containing the rotation matrix and success flag. More...
 
class  ScaledPolynomialTransform
 A 2-d coordinate transform represented by a lazy composition of an AffineTransform, a PolynomialTransform, and another AffineTransform. More...
 
class  ScaledPolynomialTransformFitter
 A fitter class for scaled polynomial transforms. More...
 
struct  ShiftRotMatrixResult
 Result of create_shift_rot_matrix() containing the rotation matrix and rotation angle. More...
 
class  SipForwardTransform
 A transform that maps pixel coordinates to intermediate world coordinates according to the SIP convention. More...
 
class  SipReverseTransform
 A transform that maps intermediate world coordinates to pixel coordinates according to the SIP convention. More...
 
class  SipTransformBase
 Base class for SIP transform objects. More...
 
struct  SortedArrayResult
 Result of create_sorted_arrays(), containing the sorted distances and array indexes. More...
 

Typedefs

typedef std::vector< RecordProxyProxyVector
 

Functions

template<typename MatchT >
afw::math::Statistics makeMatchStatistics (std::vector< MatchT > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl=afw::math::StatisticsControl())
 Compute statistics of the distance field of a match list.
 
template<typename MatchT >
afw::math::Statistics makeMatchStatisticsInPixels (afw::geom::SkyWcs const &wcs, std::vector< MatchT > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl=afw::math::StatisticsControl())
 Compute statistics of on-detector radial separation for a match list, in pixels.
 
template<typename MatchT >
afw::math::Statistics makeMatchStatisticsInRadians (afw::geom::SkyWcs const &wcs, std::vector< MatchT > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl=afw::math::StatisticsControl())
 Compute statistics of on-sky radial separation for a match list, in radians.
 
ProxyVector makeProxies (afw::table::SourceCatalog const &sourceCat, afw::geom::SkyWcs const &distortedWcs, afw::geom::SkyWcs const &tanWcs)
 
ProxyVector makeProxies (afw::table::SimpleCatalog const &posRefCat, afw::geom::SkyWcs const &tanWcs)
 
afw::table::ReferenceMatchVector matchOptimisticB (afw::table::SimpleCatalog const &posRefCat, afw::table::SourceCatalog const &sourceCat, MatchOptimisticBControl const &control, afw::geom::SkyWcs const &wcs, int posRefBegInd=0, bool verbose=false)
 Match sources to stars in a position reference catalog using optimistic pattern matching B.
 
PatternResult construct_pattern_and_shift_rot_matrix (ndarray::Array< double, 2, 1 > src_pattern_array, ndarray::Array< double, 2, 1 > src_delta_array, ndarray::Array< double, 1, 1 > src_dist_array, ndarray::Array< float, 1, 1 > dist_array, ndarray::Array< uint16_t, 2, 1 > id_array, ndarray::Array< double, 2, 1 > reference_array, size_t n_match, double max_cos_theta_shift, double max_cos_rot_sq, double max_dist_rad)
 Test an input source pattern against the reference catalog.
 
RotationTestResult test_rotation (ndarray::Array< double, 1, 1 > const &src_center, ndarray::Array< double, 1, 1 > const &ref_center, ndarray::Array< double, 1, 1 > const &src_delta, ndarray::Array< double, 1, 1 > const &ref_delta, double cos_shift, double max_cos_rot_sq)
 Test if the rotation implied between the source pattern and reference pattern is within tolerance.
 
SortedArrayResult create_sorted_arrays (ndarray::Array< double, 1, 1 > const &ref_center, ndarray::Array< double, 2, 1 > const &reference_array)
 Create arrays sorted on the distance between the center of this candidate reference object and all the reference objects.
 
std::vector< std::pair< size_t, size_t > > create_pattern_spokes (ndarray::Array< double, 1, 1 > const &src_ctr, ndarray::Array< double, 2, 1 > const &src_delta_array, ndarray::Array< double, 1, 1 > const &src_dist_array, ndarray::Array< double, 1, 1 > const &ref_ctr, Eigen::Vector3d const &proj_ref_ctr_delta, std::vector< float > const &ref_dist_array, std::vector< uint16_t > const &ref_id_array, ndarray::Array< double, 2, 1 > const &reference_array, double max_dist_rad, size_t n_match)
 Create the individual spokes that make up the pattern now that the shift and rotation are within tolerance.
 
ShiftRotMatrixResult create_shift_rot_matrix (double cos_rot_sq, Eigen::Matrix3d const &shift_matrix, ndarray::Array< double, 1, 1 > const &src_delta, ndarray::Array< double, 1, 1 > const &ref_ctr, ndarray::Array< double, 1, 1 > const &ref_delta)
 Create the complete spherical rotation matrix.
 
Eigen::Matrix3d create_spherical_rotation_matrix (Eigen::Vector3d const &rot_axis, double cos_rotation, double sin_rotion)
 Construct a generalized 3D rotation matrix about a given axis.
 
bool intermediate_verify_comparison (std::vector< Eigen::Vector3d > const &src_pattern, std::vector< Eigen::Vector3d > const &ref_pattern, Eigen::Matrix3d const &shift_rot_matrix, double max_dist_rad)
 Test the input rotation matrix by comparing the rotated src pattern to the ref pattern.
 
int check_spoke (double cos_theta_src, double sin_theta_src, ndarray::Array< double, 1, 1 > const &ref_ctr, Eigen::Vector3d const &proj_ref_ctr_delta, double proj_ref_ctr_dist_sq, std::pair< size_t, size_t > const &candidate_range, std::vector< uint16_t > 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 reference object.
 
PolynomialTransform compose (geom::AffineTransform const &t1, PolynomialTransform const &t2)
 Return a PolynomialTransform that is equivalent to the composition t1(t2())
 
PolynomialTransform compose (PolynomialTransform const &t1, geom::AffineTransform const &t2)
 Return a PolynomialTransform that is equivalent to the composition t1(t2())
 
std::shared_ptr< afw::geom::SkyWcsmakeWcs (SipForwardTransform const &sipForward, SipReverseTransform const &sipReverse, geom::SpherePoint const &skyOrigin)
 Create a new TAN SIP Wcs from a pair of SIP transforms and the sky origin.
 
std::shared_ptr< afw::geom::SkyWcstransformWcsPixels (afw::geom::SkyWcs const &wcs, geom::AffineTransform const &s)
 Create a new SkyWcs whose pixel coordinate system has been transformed via an affine transform.
 
std::shared_ptr< afw::geom::SkyWcsrotateWcsPixelsBy90 (afw::geom::SkyWcs const &wcs, int nQuarter, geom::Extent2I const &dimensions)
 Return a new SkyWcs that represents a rotation of the image it corresponds to about the image's center.
 
void wrapPolynomialTransform (WrapperCollection &wrappers)
 
void wrapScaledPolynomialTransformFitter (WrapperCollection &wrappers)
 
void wrapSipTransform (WrapperCollection &wrappers)
 
void wrapMatchOptimisticB (WrapperCollection &wrappers)
 
void wrapMakeMatchStatistics (WrapperCollection &wrappers)
 
void wrapPessimisticPatternMatcherUtils (WrapperCollection &wrappers)
 
 PYBIND11_MODULE (_measAstromLib, mod)
 
template afw::math::Statistics makeMatchStatistics< afw::table::ReferenceMatch > (std::vector< afw::table::ReferenceMatch > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl)
 
template afw::math::Statistics makeMatchStatisticsInPixels< afw::table::ReferenceMatch > (afw::geom::SkyWcs const &wcs, std::vector< afw::table::ReferenceMatch > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl)
 
template afw::math::Statistics makeMatchStatisticsInRadians< afw::table::ReferenceMatch > (afw::geom::SkyWcs const &wcs, std::vector< afw::table::ReferenceMatch > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl)
 
template afw::math::Statistics makeMatchStatistics< afw::table::SourceMatch > (std::vector< afw::table::SourceMatch > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl)
 
template afw::math::Statistics makeMatchStatisticsInPixels< afw::table::SourceMatch > (afw::geom::SkyWcs const &wcs, std::vector< afw::table::SourceMatch > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl)
 
template afw::math::Statistics makeMatchStatisticsInRadians< afw::table::SourceMatch > (afw::geom::SkyWcs const &wcs, std::vector< afw::table::SourceMatch > const &matchList, int const flags, afw::math::StatisticsControl const &sctrl)
 
std::pair< size_t, size_t > find_candidate_reference_pair_range (float src_dist, ndarray::Array< float, 1, 1 > const &ref_dist_array, double max_dist_rad)
 Find the range of reference pairs within the distance tolerance of our source pair spoke.
 
std::pair< size_t, size_t > find_candidate_reference_pair_range (float src_dist, std::vector< float > const &ref_dist_array, double max_dist_rad)
 Find the range of reference pairs within the distance tolerance of our source pair spoke.
 

Typedef Documentation

◆ ProxyVector

Definition at line 52 of file matchOptimisticB.h.

Function Documentation

◆ check_spoke()

int lsst::meas::astrom::check_spoke ( double cos_theta_src,
double sin_theta_src,
ndarray::Array< double, 1, 1 > const & ref_ctr,
Eigen::Vector3d const & proj_ref_ctr_delta,
double proj_ref_ctr_dist_sq,
std::pair< size_t, size_t > const & candidate_range,
std::vector< uint16_t > 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 reference object.

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

Parameters
[in]cos_theta_srcCosine of the angle between the current candidate source spoke and the first spoke.
[in]sin_theta_srcSine of the angle between the current candidate source spoke and the first spoke.
[in]ref_ctr3 vector of the candidate reference center
[in]proj_ref_ctr_deltaPlane projected first spoke in the reference pattern using the pattern center as normal.
[in]proj_ref_ctr_dist_sqSquared length of the projected vector.
[in]candidate_rangeMin and max index locations in ref_id_array that have pair lengths within the tolerance range.
[in]ref_id_arrayArray of id lookups into the master reference array that our center id object is paired with. uint16 type is to limit the amount of memory used and is set in the pessimistic_pattern_matcher_3d python class with reference catalogs trimmed by the matchPessimisticB runner class.
[in]reference_arrayArray of three vectors representing the locations of all reference objects.
[in]src_sin_tolSine of tolerance allowed between source and reference spoke opening angles.
Returns
ID of the candidate reference object successfully matched or -1 if no match is found.

Definition at line 368 of file pessimisticPatternMatcherUtils.cc.

371 {
372 // Loop over our candidate reference objects. candidate_range is the min
373 // and max of for pair candidates and are view into ref_id_array. Here we
374 // start from the midpoint of min and max values and step outward.
375 size_t midpoint = (candidate_range.first + candidate_range.second) / 2;
376 for (size_t idx = 0; idx < candidate_range.second - candidate_range.first; idx++) {
377 // TODO DM-33514: cleanup this loop to use an iterator that handles the "inside-out" iteration.
378 if (idx % 2 == 0) {
379 midpoint = midpoint + idx;
380 } else {
381 midpoint = midpoint - idx;
382 }
383 // Compute the delta vector from the pattern center.
384 uint16_t ref_id = ref_id_array[midpoint];
385 ndarray::Array<double, 1, 1> ref_delta = copy(reference_array[ref_id] - ref_ctr);
386
387 double ref_dot = ndarray::asEigenMatrix(ref_delta).dot(ndarray::asEigenMatrix(ref_ctr));
388 ndarray::Array<double, 1, 1> proj_ref_delta = copy(ref_delta - ref_dot * ref_ctr);
389 // Compute the cos between our "center" reference vector and the
390 // current reference candidate.
391 auto proj_ref_delta_eigen = ndarray::asEigenMatrix(proj_ref_delta);
392 double proj_delta_dist_sq = proj_ref_delta_eigen.dot(proj_ref_delta_eigen);
393 double geom_dist_ref = sqrt(proj_ref_ctr_dist_sq * proj_delta_dist_sq);
394 double cos_theta_ref = proj_ref_delta_eigen.dot(proj_ref_ctr_delta) / geom_dist_ref;
395
396 // Make sure we can safely make the comparison in case
397 // our "center" and candidate vectors are mostly aligned.
398 double cos_sq_comparison;
399 if (cos_theta_ref * cos_theta_ref < (1 - src_sin_tol * src_sin_tol)) {
400 cos_sq_comparison = (cos_theta_src - cos_theta_ref) * (cos_theta_src - cos_theta_ref) /
401 (1 - cos_theta_ref * cos_theta_ref);
402 } else {
403 cos_sq_comparison = (cos_theta_src - cos_theta_ref) * (cos_theta_src - cos_theta_ref) /
404 (src_sin_tol * src_sin_tol);
405 }
406 // Test the difference of the cosine of the reference angle against
407 // the source angle. Assumes that the delta between the two is
408 // small.
409 if (cos_sq_comparison > src_sin_tol * src_sin_tol) {
410 continue;
411 }
412 // The cosine tests the magnitude of the angle but not
413 // its direction. To do that we need to know the sine as well.
414 // This cross product calculation does that.
415 Eigen::Vector3d cross_ref = proj_ref_delta_eigen.head<3>().cross(proj_ref_ctr_delta) / geom_dist_ref;
416 double sin_theta_ref = cross_ref.dot(ndarray::asEigenMatrix(ref_ctr));
417 // Check the value of the cos again to make sure that it is not
418 // near zero.
419 double sin_comparison;
420 if (abs(cos_theta_src) < src_sin_tol) {
421 sin_comparison = (sin_theta_src - sin_theta_ref) / src_sin_tol;
422 } else {
423 sin_comparison = (sin_theta_src - sin_theta_ref) / cos_theta_ref;
424 }
425
426 // Return the correct id of the candidate we found.
427 if (abs(sin_comparison) < src_sin_tol) {
428 return ref_id;
429 }
430 }
431 return -1;
432}
Angle abs(Angle const &a)
Definition Angle.h:113

◆ compose() [1/2]

PolynomialTransform lsst::meas::astrom::compose ( geom::AffineTransform const & t1,
PolynomialTransform const & t2 )

Return a PolynomialTransform that is equivalent to the composition t1(t2())

The returned composition would be exact in ideal arithmetic, but may suffer from significant round-off error for high-order polynomials.

Definition at line 214 of file PolynomialTransform.cc.

214 {
215 typedef geom::AffineTransform AT;
216 PolynomialTransform result(t2.getOrder());
217
218 result._xCoeffs.deep() = t2._xCoeffs * t1[AT::XX] + t2._yCoeffs * t1[AT::XY];
219 result._yCoeffs.deep() = t2._xCoeffs * t1[AT::YX] + t2._yCoeffs * t1[AT::YY];
220 result._xCoeffs(0, 0) += t1[AT::X];
221 result._yCoeffs(0, 0) += t1[AT::Y];
222 return result;
223}
py::object result
Definition _schema.cc:429
An affine coordinate transformation consisting of a linear transformation and an offset.
A 2-d coordinate transform represented by a pair of standard polynomials (one for each coordinate).

◆ compose() [2/2]

PolynomialTransform lsst::meas::astrom::compose ( PolynomialTransform const & t1,
geom::AffineTransform const & t2 )

Return a PolynomialTransform that is equivalent to the composition t1(t2())

The returned composition would be exact in ideal arithmetic, but may suffer from significant round-off error for high-order polynomials.

Definition at line 225 of file PolynomialTransform.cc.

225 {
226 typedef geom::AffineTransform AT;
227 int const order = t1.getOrder();
228 if (order < 1) {
229 PolynomialTransform t1a(1);
230 t1a._xCoeffs(0, 0) = t1._xCoeffs(0, 0);
231 t1a._yCoeffs(0, 0) = t1._yCoeffs(0, 0);
232 return compose(t1a, t2);
233 }
234 detail::BinomialMatrix binomial(order);
235 // For each of these, (e.g.) a[n] == pow(a, n)
236 auto const t2u = detail::computePowers(t2[AT::X], order);
237 auto const t2v = detail::computePowers(t2[AT::Y], order);
238 auto const t2uu = detail::computePowers(t2[AT::XX], order);
239 auto const t2uv = detail::computePowers(t2[AT::XY], order);
240 auto const t2vu = detail::computePowers(t2[AT::YX], order);
241 auto const t2vv = detail::computePowers(t2[AT::YY], order);
242 PolynomialTransform result(order);
243 for (int p = 0; p <= order; ++p) {
244 for (int m = 0; m <= p; ++m) {
245 for (int j = 0; j <= m; ++j) {
246 for (int q = 0; p + q <= order; ++q) {
247 for (int n = 0; n <= q; ++n) {
248 for (int k = 0; k <= n; ++k) {
249 double z = binomial(p, m) * t2u[p - m] * binomial(m, j) * t2uu[j] * t2uv[m - j] *
250 binomial(q, n) * t2v[q - n] * binomial(n, k) * t2vu[k] * t2vv[n - k];
251 result._xCoeffs(j + k, m + n - j - k) += t1._xCoeffs(p, q) * z;
252 result._yCoeffs(j + k, m + n - j - k) += t1._yCoeffs(p, q) * z;
253 } // k
254 } // n
255 } // q
256 } // j
257 } // m
258 } // p
259 return result;
260}
double z
Definition Match.cc:44
int m
Definition SpanSet.cc:48
table::Key< int > order

◆ construct_pattern_and_shift_rot_matrix()

PatternResult lsst::meas::astrom::construct_pattern_and_shift_rot_matrix ( ndarray::Array< double, 2, 1 > src_pattern_array,
ndarray::Array< double, 2, 1 > src_delta_array,
ndarray::Array< double, 1, 1 > src_dist_array,
ndarray::Array< float, 1, 1 > dist_array,
ndarray::Array< uint16_t, 2, 1 > id_array,
ndarray::Array< double, 2, 1 > reference_array,
size_t n_match,
double max_cos_theta_shift,
double max_cos_rot_sq,
double max_dist_rad )

Test an input source pattern against the reference catalog.

Parameters
[in]src_pattern_arraySub selection of source 3 vectors to create a pattern from.
[in]src_delta_arrayDeltas of pairs from the source pattern center.
[in]src_dist_arrayDistances of the pairs in src_delta_array.
[in]dist_arraySet of all distances between pairs of reference objects.
[in]id_arrayPairs of indices creating the reference pairs. Sorted the same as dist_array.
[in]reference_arraySet of all reference points.
[in]n_matchNumber of points to attempt to create a pattern from. Must be <= len(src_pattern_array)
[in]max_cos_theta_shiftMaximum shift allowed between two patterns' centers. Equivalent to the shift on the sky.
[in]max_cos_rot_sqMaximum rotation between two patterns that have been shifted to have their centers on top of each other.
[in]max_dist_radMaximum difference allowed between the source and reference pair distances to consider the reference pair a candidate for the source pair. Also sets the tolerance between the opening angles of the spokes when compared to the reference.
Returns
The result of the pattern matching.

Definition at line 40 of file pessimisticPatternMatcherUtils.cc.

44 {
45 // Our first test. We search the reference dataset for pairs that have the same length as our first source
46 // pairs to within plus/minus the max_dist tolerance.
47 std::pair<size_t, size_t> candidate_range =
48 find_candidate_reference_pair_range(src_dist_array[0], dist_array, max_dist_rad);
49 size_t ref_dist_idx = (candidate_range.first + candidate_range.second) / 2;
50
51 // Start our loop over the candidate reference objects. Looping from the inside (minimum difference to our
52 // source dist) to the outside.
53 for (size_t idx = 0; idx < candidate_range.second - candidate_range.first; idx++) {
54 // TODO DM-33514: cleanup this loop to use an iterator that handles the "inside-out" iteration.
55 if (idx % 2 == 0) {
56 ref_dist_idx = ref_dist_idx + idx;
57 } else {
58 ref_dist_idx = ref_dist_idx - idx;
59 }
60 // We have two candidates for which reference object corresponds with the source at the center of our
61 // pattern. As such we loop over and test both possibilities.
62 ndarray::Array<uint16_t, 1, 1> tmp_ref_pair_list = id_array[ref_dist_idx];
63 for (uint16_t ref_pair_idx = 0; ref_pair_idx < 2; ref_pair_idx++) {
64 uint16_t ref_id = id_array[ref_dist_idx][ref_pair_idx];
66 // Test the angle between our candidate ref center and the source center of our pattern. This
67 // angular distance also defines the shift we will later use.
68 ndarray::Array<double, 1, 1> ref_center = reference_array[ref_id];
69 double cos_shift =
70 ndarray::asEigenMatrix(src_pattern_array[0]).dot(ndarray::asEigenMatrix(ref_center));
71 if (cos_shift < max_cos_theta_shift) {
72 continue;
73 }
74 // We can now append this one as a candidate.
75 candidate_pairs.push_back(std::make_pair(ref_id, 0));
76 ndarray::Array<double, 1, 1> ref_delta;
77 // Test to see which reference object to use in the pair.
78 if (ref_id == *tmp_ref_pair_list.begin()) {
79 candidate_pairs.push_back(std::make_pair(tmp_ref_pair_list[1], 1));
80 ref_delta = copy(reference_array[tmp_ref_pair_list[1]] - ref_center);
81 } else {
82 candidate_pairs.push_back(std::make_pair(tmp_ref_pair_list[0], 1));
83 ref_delta = copy(reference_array[tmp_ref_pair_list[0]] - ref_center);
84 }
85 // For dense fields it will be faster to compute the absolute rotation this pair suggests first
86 // rather than saving it after all the spokes are found. We then compute the cos^2 of the rotation
87 // and first part of the rotation matrix from source to reference frame.
88 RotationTestResult test_rot_result =
89 test_rotation(src_pattern_array[0], ref_center, src_delta_array[0], ref_delta, cos_shift,
90 max_cos_rot_sq);
91 if (!test_rot_result.success) {
92 continue;
93 }
94 // Now that we have a candidate first spoke and reference pattern center, we mask our future
95 // search to only those pairs that contain our candidate reference center.
96 SortedArrayResult sorted_array_struct = create_sorted_arrays(ref_center, reference_array);
97 // Now we feed this sub data to match the spokes of our pattern.
99 create_pattern_spokes(src_pattern_array[0], src_delta_array, src_dist_array, ref_center,
100 test_rot_result.proj_ref_ctr_delta, sorted_array_struct.dists,
101 sorted_array_struct.ids, reference_array, max_dist_rad, n_match);
102 // If we don't find enough candidates we can continue to the next reference center pair.
103 if (pattern_spokes.size() < n_match - 2) {
104 continue;
105 }
106 // If we have the right number of matched ids we store these.
107 candidate_pairs.reserve(candidate_pairs.size() + pattern_spokes.size());
108 candidate_pairs.insert(candidate_pairs.end(), pattern_spokes.begin(), pattern_spokes.end());
109 // We can now create our full matrix for both the shift and rotation. The shift aligns
110 // the pattern centers, while the rotation rotates the spokes on top of each other.
111 ShiftRotMatrixResult shift_rot_result =
112 create_shift_rot_matrix(test_rot_result.cos_rot_sq, test_rot_result.shift_matrix,
113 src_delta_array[0], ref_center, ref_delta);
114
115 // concatenate our final patterns.
116 Eigen::Matrix3d shift_rot_matrix = shift_rot_result.shift_rot_matrix;
119 ref_pattern.reserve(n_match);
120 src_pattern.reserve(n_match);
121 for (size_t idx = 0; idx < n_match; idx++) {
122 ref_pattern[idx] = ndarray::asEigenMatrix(reference_array[candidate_pairs[idx].first]);
123 src_pattern[idx] = ndarray::asEigenMatrix(src_pattern_array[candidate_pairs[idx].second]);
124 }
125 // TODO: DM-32985 Implement least squares linear fitting here. Look at python example linked in
126 // ticket for how to implement.
127 // Test that the all points in each pattern are within tolerance after shift/rotation.
128 bool passed =
129 intermediate_verify_comparison(src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad);
130 if (passed) {
131 return PatternResult(candidate_pairs, shift_rot_matrix, cos_shift, shift_rot_result.sin_rot);
132 }
133 }
134 }
135 // failed fit
136 return PatternResult();
137}
T begin(T... args)
T end(T... args)
T insert(T... args)
T make_pair(T... args)
RotationTestResult test_rotation(ndarray::Array< double, 1, 1 > const &src_center, ndarray::Array< double, 1, 1 > const &ref_center, ndarray::Array< double, 1, 1 > const &src_delta, ndarray::Array< double, 1, 1 > const &ref_delta, double cos_shift, double max_cos_rot_sq)
Test if the rotation implied between the source pattern and reference pattern is within tolerance.
bool intermediate_verify_comparison(std::vector< Eigen::Vector3d > const &src_pattern, std::vector< Eigen::Vector3d > const &ref_pattern, Eigen::Matrix3d const &shift_rot_matrix, double max_dist_rad)
Test the input rotation matrix by comparing the rotated src pattern to the ref pattern.
SortedArrayResult create_sorted_arrays(ndarray::Array< double, 1, 1 > const &ref_center, ndarray::Array< double, 2, 1 > const &reference_array)
Create arrays sorted on the distance between the center of this candidate reference object and all th...
std::pair< size_t, size_t > find_candidate_reference_pair_range(float src_dist, ndarray::Array< float, 1, 1 > const &ref_dist_array, double max_dist_rad)
Find the range of reference pairs within the distance tolerance of our source pair spoke.
ShiftRotMatrixResult create_shift_rot_matrix(double cos_rot_sq, Eigen::Matrix3d const &shift_matrix, ndarray::Array< double, 1, 1 > const &src_delta, ndarray::Array< double, 1, 1 > const &ref_ctr, ndarray::Array< double, 1, 1 > const &ref_delta)
Create the complete spherical rotation matrix.
std::vector< std::pair< size_t, size_t > > create_pattern_spokes(ndarray::Array< double, 1, 1 > const &src_ctr, ndarray::Array< double, 2, 1 > const &src_delta_array, ndarray::Array< double, 1, 1 > const &src_dist_array, ndarray::Array< double, 1, 1 > const &ref_ctr, Eigen::Vector3d const &proj_ref_ctr_delta, std::vector< float > const &ref_dist_array, std::vector< uint16_t > const &ref_id_array, ndarray::Array< double, 2, 1 > const &reference_array, double max_dist_rad, size_t n_match)
Create the individual spokes that make up the pattern now that the shift and rotation are within tole...
T push_back(T... args)
T reserve(T... args)
T size(T... args)

◆ create_pattern_spokes()

std::vector< std::pair< size_t, size_t > > lsst::meas::astrom::create_pattern_spokes ( ndarray::Array< double, 1, 1 > const & src_ctr,
ndarray::Array< double, 2, 1 > const & src_delta_array,
ndarray::Array< double, 1, 1 > const & src_dist_array,
ndarray::Array< double, 1, 1 > const & ref_ctr,
Eigen::Vector3d const & proj_ref_ctr_delta,
std::vector< float > const & ref_dist_array,
std::vector< uint16_t > const & ref_id_array,
ndarray::Array< double, 2, 1 > const & reference_array,
double max_dist_rad,
size_t n_match )

Create the individual spokes that make up the pattern now that the shift and rotation are within tolerance.

If we can't create a valid pattern we exit early with a partial result.

Parameters
[in]src_ctr3 vector of the source pinwheel center
[in]src_delta_arrayArray of 3 vector deltas between the source center and the pairs that make up the remaining spokes of the pinwheel.
[in]src_dist_arrayArray of the distances of each src_delta in the pinwheel.
[in]ref_ctr3 vector of the candidate reference center.
[in]proj_ref_ctr_deltaPlane projected 3 vector formed from the center point of the candidate pin-wheel and the second point in the pattern to create the first spoke pair. This is the candidate pair that was matched in the main _construct_pattern_and_shift_rot_matrix loop.
[in]ref_dist_arrayArray of vector distances for each of the reference pairs
[in]ref_id_arrayArray of id lookups into the master reference array that our center id object is paired with.
[in]reference_arrayFull 3 vector data for the reference catalog.
[in]max_dist_radMaximum search radius for distances.
[in]n_matchNumber of source deltas that must be matched into the reference deltas in order to consider this a successful pattern match.
Returns
Return pairs of reference ids and their matched src ids.

Append the successful indices to our list. The src_idx needs an extra iteration to skip the first and second source objects.

Append the successful indices to our list. The src_idx needs an extra iteration to skip the first and second source objects.

Definition at line 278 of file pessimisticPatternMatcherUtils.cc.

283 {
284 // Struct where we will be putting our results.
286
287 // Counter for number of spokes we failed to find a reference
288 // candidate for. We break the loop if we haven't found enough.
289 size_t n_fail = 0;
290
291 // Plane project the center/first spoke of the source pattern using
292 // the center vector of the pattern as normal.
293 auto src_ctr_eigen = ndarray::asEigenMatrix(src_ctr);
294 double src_delta_ctr_dot = ndarray::asEigenMatrix(src_delta_array[0]).dot(src_ctr_eigen);
295
296 ndarray::Array<double, 1, 1> proj_src_ctr_delta = copy(src_delta_array[0] - src_delta_ctr_dot * src_ctr);
297 auto proj_src_ctr_delta_eigen = ndarray::asEigenMatrix(proj_src_ctr_delta);
298 double proj_src_ctr_dist_sq = proj_src_ctr_delta_eigen.dot(proj_src_ctr_delta_eigen);
299
300 // Pre - compute the squared length of the projected reference vector.
301 double proj_ref_ctr_dist_sq = proj_ref_ctr_delta.dot(proj_ref_ctr_delta);
302
303 // Value of sin where sin(theta) ~= theta to within 0.1%. Used to make
304 // sure we are still in small angle approximation.
305 double max_sin_tol = 0.0447;
306 // Loop over the source pairs.
307 for (size_t src_idx = 1; src_idx < src_dist_array.size(); src_idx++) {
308 if (n_fail > src_dist_array.size() - (n_match - 1)) {
309 break;
310 }
311 // Find the reference pairs that include our candidate pattern center
312 // and sort them in increasing delta. Check this first so we don't
313 // compute anything else if no candidates exist.
314 std::pair<size_t, size_t> candidate_range =
315 find_candidate_reference_pair_range(src_dist_array[src_idx], ref_dist_array, max_dist_rad);
316 if (candidate_range.first == candidate_range.second) {
317 n_fail++;
318 continue;
319 }
320
321 // Given our length tolerance we can use it to compute a tolerance
322 // on the angle between our spoke.
323 double src_sin_tol = max_dist_rad / (src_dist_array[src_idx] + max_dist_rad);
324
325 // Test if the small angle approximation will still hold. This is
326 // defined as when sin(theta) ~= theta to within 0.1% of each
327 // other. If the implied opening angle is too large we set it to
328 // the 0.1% threshold.
329 if (src_sin_tol > max_sin_tol) {
330 src_sin_tol = max_sin_tol;
331 }
332
333 // Plane project the candidate source spoke and compute the cosine
334 // and sine of the opening angle.
335 double proj_src_delta_dot = ndarray::asEigenMatrix(src_delta_array[src_idx]).dot(src_ctr_eigen);
336
337 ndarray::Array<double, 1, 1> proj_src_delta =
338 copy(src_delta_array[src_idx] - proj_src_delta_dot * src_ctr);
339 auto proj_src_delta_eigen = ndarray::asEigenMatrix(proj_src_delta);
340 double geom_dist_src = sqrt(proj_src_delta_eigen.dot(proj_src_delta_eigen) * proj_src_ctr_dist_sq);
341
342 // Compute cosine and sine of the delta vector opening angle.
343 double cos_theta_src = proj_src_delta_eigen.dot(proj_src_ctr_delta_eigen) / geom_dist_src;
344 Eigen::Vector3d cross_src =
345 proj_src_delta_eigen.head<3>().cross(proj_src_ctr_delta_eigen.head<3>()) / geom_dist_src;
346 double sin_theta_src = cross_src.dot(src_ctr_eigen);
347
348 // Test the spokes and return the id of the reference object.
349 // Return -1 if no match is found.
350 int ref_id =
351 check_spoke(cos_theta_src, sin_theta_src, ref_ctr, proj_ref_ctr_delta, proj_ref_ctr_dist_sq,
352 candidate_range, ref_id_array, reference_array, src_sin_tol);
353 if (ref_id < 0) {
354 n_fail++;
355 continue;
356 }
357
360 output_spokes.push_back(std::make_pair(static_cast<size_t>(ref_id), src_idx + 1));
361 if (output_spokes.size() >= n_match - 2) {
362 break;
363 }
364 }
365 return output_spokes;
366}

◆ create_shift_rot_matrix()

ShiftRotMatrixResult lsst::meas::astrom::create_shift_rot_matrix ( double cos_rot_sq,
Eigen::Matrix3d const & shift_matrix,
ndarray::Array< double, 1, 1 > const & src_delta,
ndarray::Array< double, 1, 1 > const & ref_ctr,
ndarray::Array< double, 1, 1 > const & ref_delta )

Create the complete spherical rotation matrix.

Parameters
[in]cos_rot_sqcosine of the rotation needed to align our source and reference candidate patterns.
[in]shift_matrix3x3 rotation matrix for shifting the source pattern center on top of the candidate reference pattern center.
[in]src_delta3 vector delta representing the first spoke of the source pattern
[in]ref_ctr3 vector on the unit-sphere representing the center of our reference pattern.
[in]ref_delta3 vector delta made by the first pair of the reference pattern.
Returns
Struct containing constructed matrix and implied rotation.

Definition at line 237 of file pessimisticPatternMatcherUtils.cc.

240 {
241 double cos_rot = sqrt(cos_rot_sq);
242 Eigen::Vector3d src_delta_eigen = ndarray::asEigenMatrix(src_delta);
243 Eigen::Vector3d rot_src_delta = shift_matrix * src_delta_eigen;
244 Eigen::Vector3d tmp_cross = rot_src_delta.cross(ndarray::asEigenMatrix(ref_delta).head<3>());
245 double delta_dot_cross = tmp_cross.dot(ndarray::asEigenMatrix(ref_ctr));
246
247 double sin_rot = sgn(delta_dot_cross) * sqrt(1 - cos_rot_sq);
248 Eigen::Matrix3d rot_matrix =
249 create_spherical_rotation_matrix(ndarray::asEigenMatrix(ref_ctr), cos_rot, sin_rot);
250
251 Eigen::Matrix3d shift_rot_matrix = rot_matrix * shift_matrix;
252
254 result.sin_rot = sin_rot;
255 result.shift_rot_matrix = shift_rot_matrix;
256 return result;
257}
Eigen::Matrix3d create_spherical_rotation_matrix(Eigen::Vector3d const &rot_axis, double cos_rotation, double sin_rotion)
Construct a generalized 3D rotation matrix about a given axis.
Result of create_shift_rot_matrix() containing the rotation matrix and rotation angle.
double sin_rot
Rotation that makes up the matrix.

◆ create_sorted_arrays()

SortedArrayResult lsst::meas::astrom::create_sorted_arrays ( ndarray::Array< double, 1, 1 > const & ref_center,
ndarray::Array< double, 2, 1 > const & reference_array )

Create arrays sorted on the distance between the center of this candidate reference object and all the reference objects.

Parameters
[in]ref_centerCenter point of the candidate pattern.
[in]reference_arrayArray of all reference object points.
Returns
Sorted distances and indexes of reference objects.

Definition at line 139 of file pessimisticPatternMatcherUtils.cc.

140 {
142 // NOTE: this algorithm is quadratic in the length of reference_array. It might be worth using std::sort
143 // instead of this approach, if reference_array is long, but the length at which that matters should
144 // be tested because it would require caching of distances.
145 for (uint16_t idx = 0; idx < reference_array.getShape()[0]; idx++) {
146 Eigen::Vector3d diff = ndarray::asEigenMatrix(copy(reference_array[idx] - ref_center));
147 double dist = sqrt(diff.dot(diff));
148 auto dists_itr = std::lower_bound(result.dists.begin(), result.dists.end(), dist);
149 auto ids_itr = result.ids.begin() + (dists_itr - result.dists.begin());
150 result.ids.insert(ids_itr, idx);
151 result.dists.insert(dists_itr, dist);
152 }
153 return result;
154}
T lower_bound(T... args)
Result of create_sorted_arrays(), containing the sorted distances and array indexes.

◆ create_spherical_rotation_matrix()

Eigen::Matrix3d lsst::meas::astrom::create_spherical_rotation_matrix ( Eigen::Vector3d const & rot_axis,
double cos_rotation,
double sin_rotion )

Construct a generalized 3D rotation matrix about a given axis.

param[in] rot_axis 3 vector defining the axis to rotate about. param[in] cos_rotation cosine of the rotation angle. param[in] sin_rotion sine of the rotation angle.

Returns
3x3 spherical, rotation matrix.

Definition at line 223 of file pessimisticPatternMatcherUtils.cc.

224 {
225 Eigen::Matrix3d rot_cross_matrix;
226 // clang-format off
227 rot_cross_matrix << 0., -rot_axis[2], rot_axis[1],
228 rot_axis[2], 0., -rot_axis[0],
229 -rot_axis[1], rot_axis[0], 0.;
230 // clang-format on
231 Eigen::Matrix3d shift_matrix = cos_rotation * Eigen::Matrix3d::Identity() +
232 sin_rotation * rot_cross_matrix +
233 (1 - cos_rotation) * rot_axis * rot_axis.transpose();
234 return shift_matrix;
235}

◆ find_candidate_reference_pair_range() [1/2]

std::pair< size_t, size_t > lsst::meas::astrom::find_candidate_reference_pair_range ( float src_dist,
ndarray::Array< float, 1, 1 > const & ref_dist_array,
double max_dist_rad )

Find the range of reference pairs within the distance tolerance of our source pair spoke.

Parameters
[in]src_distValue of the distance we would like to search for in the reference array in radians.
[in]ref_dist_arraysorted array of distances in radians.
[in]max_dist_radmaximum plus/minus search to find in the reference array in radians.
Returns
pair of indices for the min and max range of indices spanning src_dist +/- max_dist_rad to search.

Definition at line 156 of file pessimisticPatternMatcherUtils.cc.

157 {
158 auto itr =
159 std::lower_bound(ref_dist_array.begin(), ref_dist_array.end(), src_dist - max_dist_rad - 1e-16);
160 auto itrEnd =
161 std::upper_bound(ref_dist_array.begin(), ref_dist_array.end(), src_dist + max_dist_rad + 1e-16);
162
163 size_t startIdx = itr - ref_dist_array.begin();
164 size_t endIdx = itrEnd - ref_dist_array.begin();
165 return std::make_pair(startIdx, endIdx);
166}
T upper_bound(T... args)

◆ find_candidate_reference_pair_range() [2/2]

std::pair< size_t, size_t > lsst::meas::astrom::find_candidate_reference_pair_range ( float src_dist,
std::vector< float > const & ref_dist_array,
double max_dist_rad )

Find the range of reference pairs within the distance tolerance of our source pair spoke.

Parameters
[in]src_distValue of the distance we would like to search for in the reference array in radians.
[in]ref_dist_arraysorted array of distances in radians.
[in]max_dist_radmaximum plus/minus search to find in the reference array in radians.
Returns
pair of indices for the min and max range of indices spanning src_dist +/- max_dist_rad to search.

Definition at line 168 of file pessimisticPatternMatcherUtils.cc.

170 {
171 auto itr =
172 std::lower_bound(ref_dist_array.begin(), ref_dist_array.end(), src_dist - max_dist_rad - 1e-16);
173 auto itrEnd =
174 std::upper_bound(ref_dist_array.begin(), ref_dist_array.end(), src_dist + max_dist_rad + 1e-16);
175
176 size_t startIdx = itr - ref_dist_array.begin();
177 size_t endIdx = itrEnd - ref_dist_array.begin();
178 return std::make_pair(startIdx, endIdx);
179}

◆ intermediate_verify_comparison()

bool lsst::meas::astrom::intermediate_verify_comparison ( std::vector< Eigen::Vector3d > const & src_pattern,
std::vector< Eigen::Vector3d > const & ref_pattern,
Eigen::Matrix3d const & shift_rot_matrix,
double max_dist_rad )

Test the input rotation matrix by comparing the rotated src pattern to the ref pattern.

If every point in the pattern after rotation is within a distance of max_dist_rad to its candidate point in the other pattern, we return True.

Parameters
[in]src_patternVector of 3 vectors representing the points that make up our source pinwheel pattern.
[in]ref_patternVector of 3 vectors representing our candidate reference pinwheel pattern.
[in]shift_rot_matrix3x3 rotation matrix that takes the source objects and rotates them onto the frame of the reference objects.
[in]max_dist_radMaximum distance allowed to consider two objects the same.
Returns
Whether or not the rotation matrix gives results within tolerances.

Definition at line 259 of file pessimisticPatternMatcherUtils.cc.

261 {
262 double max_dist_sq = max_dist_rad * max_dist_rad;
263 bool passed = true;
264 auto iSrc = src_pattern.begin();
265 auto iRef = ref_pattern.begin();
266 for (auto iSrc = src_pattern.begin(), iRef = ref_pattern.begin();
267 iSrc != src_pattern.end() && iRef != ref_pattern.end(); iSrc++, iRef++) {
268 Eigen::Vector3d rot_src_vect = shift_rot_matrix * *iSrc;
269 Eigen::Vector3d diff_vect = rot_src_vect - *iRef;
270 if (max_dist_sq < diff_vect.dot(diff_vect)) {
271 passed = false;
272 break;
273 }
274 }
275 return passed;
276}

◆ makeMatchStatistics()

template<typename MatchT >
afw::math::Statistics lsst::meas::astrom::makeMatchStatistics ( std::vector< MatchT > const & matchList,
int const flags,
afw::math::StatisticsControl const & sctrl = afw::math::StatisticsControl() )

Compute statistics of the distance field of a match list.

Parameters
[in]matchListlist of matchList between reference objects and sources; fields read:
  • distance: distance between source and reference object, in arbitrary units; the resulting statistics have the same units as distance
[in]flagswhat to calculate; OR constants such as lsst::afw::math::MEAN, MEANCLIP, STDDEV, MEDIAN, defined in lsst/afw/math/Statitics.h's Property enum
[in]sctrlstatistics configuration

Definition at line 33 of file makeMatchStatistics.cc.

34 {
35 if (matchList.empty()) {
36 throw LSST_EXCEPT(pexExcept::RuntimeError, "matchList is empty");
37 }
39 val.reserve(matchList.size());
40
41 for (auto const& match : matchList) {
42 val.push_back(match.distance);
43 }
44 return afw::math::makeStatistics(val, flags, sctrl);
45}
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition Exception.h:48
Reports errors that are due to events beyond the control of the program.
Definition Runtime.h:104
T empty(T... args)
ImageT val
Definition CR.cc:146

◆ makeMatchStatistics< afw::table::ReferenceMatch >()

template afw::math::Statistics lsst::meas::astrom::makeMatchStatistics< afw::table::ReferenceMatch > ( std::vector< afw::table::ReferenceMatch > const & matchList,
int const flags,
afw::math::StatisticsControl const & sctrl )

◆ makeMatchStatistics< afw::table::SourceMatch >()

template afw::math::Statistics lsst::meas::astrom::makeMatchStatistics< afw::table::SourceMatch > ( std::vector< afw::table::SourceMatch > const & matchList,
int const flags,
afw::math::StatisticsControl const & sctrl )

◆ makeMatchStatisticsInPixels()

template<typename MatchT >
afw::math::Statistics lsst::meas::astrom::makeMatchStatisticsInPixels ( afw::geom::SkyWcs const & wcs,
std::vector< MatchT > const & matchList,
int const flags,
afw::math::StatisticsControl const & sctrl = afw::math::StatisticsControl() )

Compute statistics of on-detector radial separation for a match list, in pixels.

Parameters
[in]wcsWCS describing pixel to sky transformation
[in]matchListlist of matchList between reference objects and sources; fields read:
  • first: reference object; only the coord is read
  • second: source; only the centroid is read
[in]flagswhat to calculate; OR constants such as lsst::afw::math::MEAN, MEANCLIP, STDDEV, MEDIAN, defined in lsst/afw/math/Statitics.h's Property enum
[in]sctrlstatistics configuration

Definition at line 48 of file makeMatchStatistics.cc.

50 {
51 if (matchList.empty()) {
52 throw LSST_EXCEPT(pexExcept::RuntimeError, "matchList is empty");
53 }
55 val.reserve(matchList.size());
56
57 for (auto const& match : matchList) {
58 auto refPtr = match.first;
59 auto srcPtr = match.second;
60 auto srcX = srcPtr->getX();
61 auto srcY = srcPtr->getY();
62 auto refPos = wcs.skyToPixel(refPtr->getCoord());
63 auto refX = refPos[0];
64 auto refY = refPos[1];
65 val.push_back(::hypot(srcX - refX, srcY - refY));
66 }
67 return afw::math::makeStatistics(val, flags, sctrl);
68}

◆ makeMatchStatisticsInPixels< afw::table::ReferenceMatch >()

◆ makeMatchStatisticsInPixels< afw::table::SourceMatch >()

◆ makeMatchStatisticsInRadians()

template<typename MatchT >
afw::math::Statistics lsst::meas::astrom::makeMatchStatisticsInRadians ( afw::geom::SkyWcs const & wcs,
std::vector< MatchT > const & matchList,
int const flags,
afw::math::StatisticsControl const & sctrl = afw::math::StatisticsControl() )

Compute statistics of on-sky radial separation for a match list, in radians.

Parameters
[in]wcsWCS describing pixel to sky transformation
[in]matchListlist of matchList between reference objects and sources; fields read:
  • first: reference object; only the coord is read
  • second: source; only the centroid is read
[in]flagswhat to calculate; OR constants such as lsst::afw::math::MEAN, MEANCLIP, STDDEV, MEDIAN, defined in lsst/afw/math/Statitics.h's Property enum
[in]sctrlstatistics configuration

Definition at line 71 of file makeMatchStatistics.cc.

73 {
74 if (matchList.empty()) {
75 throw LSST_EXCEPT(pexExcept::RuntimeError, "matchList is empty");
76 }
78 val.reserve(matchList.size());
79
80 for (auto const& match : matchList) {
81 auto refPtr = match.first;
82 auto srcPtr = match.second;
83 auto refCoord = refPtr->getCoord();
84 auto srcCoord = wcs.pixelToSky(srcPtr->getCentroid());
85 auto angSep = refCoord.separation(srcCoord);
86 val.push_back(angSep.asRadians());
87 }
88 return afw::math::makeStatistics(val, flags, sctrl);
89}

◆ makeMatchStatisticsInRadians< afw::table::ReferenceMatch >()

template afw::math::Statistics lsst::meas::astrom::makeMatchStatisticsInRadians< afw::table::ReferenceMatch > ( afw::geom::SkyWcs const & wcs,
std::vector< afw::table::ReferenceMatch > const & matchList,
int const flags,
afw::math::StatisticsControl const & sctrl )

◆ makeMatchStatisticsInRadians< afw::table::SourceMatch >()

template afw::math::Statistics lsst::meas::astrom::makeMatchStatisticsInRadians< afw::table::SourceMatch > ( afw::geom::SkyWcs const & wcs,
std::vector< afw::table::SourceMatch > const & matchList,
int const flags,
afw::math::StatisticsControl const & sctrl )

◆ makeProxies() [1/2]

ProxyVector lsst::meas::astrom::makeProxies ( afw::table::SimpleCatalog const & posRefCat,
afw::geom::SkyWcs const & tanWcs )

Definition at line 448 of file matchOptimisticB.cc.

448 {
449 auto coordKey = afwTable::CoordKey(posRefCat.getSchema()["coord"]);
450 ProxyVector r;
451 r.reserve(posRefCat.size());
452 for (afwTable::SimpleCatalog::const_iterator posRefPtr = posRefCat.begin(); posRefPtr != posRefCat.end();
453 ++posRefPtr) {
454 r.push_back(RecordProxy(posRefPtr, tanWcs.skyToPixel(posRefPtr->get(coordKey))));
455 }
456 return r;
457}
A FunctorKey used to get or set celestial coordinates from a pair of lsst::geom::Angle keys.
Definition aggregates.h:292
A wrapper around a SimpleRecord or SourceRecord that allows us to record a pixel position in a way th...

◆ makeProxies() [2/2]

ProxyVector lsst::meas::astrom::makeProxies ( afw::table::SourceCatalog const & sourceCat,
afw::geom::SkyWcs const & distortedWcs,
afw::geom::SkyWcs const & tanWcs )

Definition at line 436 of file matchOptimisticB.cc.

437 {
438 ProxyVector r;
439 r.reserve(sourceCat.size());
440 for (afwTable::SourceCatalog::const_iterator sourcePtr = sourceCat.begin(); sourcePtr != sourceCat.end();
441 ++sourcePtr) {
442 r.push_back(
443 RecordProxy(sourcePtr, tanWcs.skyToPixel(distortedWcs.pixelToSky(sourcePtr->getCentroid()))));
444 }
445 return r;
446}

◆ makeWcs()

std::shared_ptr< afw::geom::SkyWcs > lsst::meas::astrom::makeWcs ( SipForwardTransform const & sipForward,
SipReverseTransform const & sipReverse,
geom::SpherePoint const & skyOrigin )

Create a new TAN SIP Wcs from a pair of SIP transforms and the sky origin.

Parameters
[in]sipForwardMapping from pixel coordinates to intermediate world coordinates.
[in]sipReverseMapping from intermediate world coordinates to pixel coordinates.
[in]skyOriginICRS position of the gnomonic projection that maps sky coordinates to intermediate world coordinates (CRVAL).
Exceptions
pex::exceptions::InvalidParameterErrorif the forward and reverse SIP transforms have different CRPIX values or CD matrices.

Definition at line 148 of file SipTransform.cc.

150 {
151 if (!sipForward.getPixelOrigin().asEigen().isApprox(sipReverse.getPixelOrigin().asEigen())) {
153 oss << "SIP forward and reverse transforms have inconsistent CRPIX: " << sipForward.getPixelOrigin()
154 << " != " << sipReverse.getPixelOrigin();
156 }
157 if (!sipForward.getCdMatrix().getMatrix().isApprox(sipReverse.getCdMatrix().getMatrix())) {
159 oss << "SIP forward and reverse transforms have inconsistent CD matrix: " << sipForward.getCdMatrix()
160 << "\n!=\n"
161 << sipReverse.getCdMatrix();
163 }
164 Eigen::MatrixXd sipA(ndarray::asEigenMatrix(sipForward.getPoly().getXCoeffs()));
165 Eigen::MatrixXd sipB(ndarray::asEigenMatrix(sipForward.getPoly().getYCoeffs()));
166 Eigen::MatrixXd sipAP(ndarray::asEigenMatrix(sipReverse.getPoly().getXCoeffs()));
167 Eigen::MatrixXd sipBP(ndarray::asEigenMatrix(sipReverse.getPoly().getYCoeffs()));
168
169 return afw::geom::makeTanSipWcs(sipForward.getPixelOrigin(), skyOrigin,
170 sipForward.getCdMatrix().getMatrix(), sipA, sipB, sipAP, sipBP);
171}
Reports invalid arguments.
Definition Runtime.h:66
T str(T... args)

◆ matchOptimisticB()

afwTable::ReferenceMatchVector lsst::meas::astrom::matchOptimisticB ( afw::table::SimpleCatalog const & posRefCat,
afw::table::SourceCatalog const & sourceCat,
MatchOptimisticBControl const & control,
afw::geom::SkyWcs const & wcs,
int posRefBegInd = 0,
bool verbose = false )

Match sources to stars in a position reference catalog using optimistic pattern matching B.

Optimistic Pattern Matching is described in V. Tabur 2007, PASA, 24, 189-198 "Fast Algorithms for Matching CCD Images to a Stellar Catalogue"

Parameters
[in]posRefCatcatalog of position reference stars; fields read:
  • "coord"
  • control.refFluxField
[in]sourceCatcatalog of detected sources; fields read:
  • "Centroid_x"
  • "Centroid_y"
  • control.refFluxField
[in]wcsestimated WCS
[in]controlcontrol object
[in]posRefBegIndindex of first start to use in posRefCat
[in]verbosetrue to print diagnostic information to std::cout
Returns
a list of matches; the d field may be set, but you should not rely on it

Definition at line 459 of file matchOptimisticB.cc.

463 {
464 control.validate();
465 if (posRefCat.empty()) {
466 throw LSST_EXCEPT(pexExcept::InvalidParameterError, "no entries in posRefCat");
467 }
468 if (sourceCat.empty()) {
469 throw LSST_EXCEPT(pexExcept::InvalidParameterError, "no entries in sourceCat");
470 }
471 if (posRefBegInd < 0) {
472 throw LSST_EXCEPT(pexExcept::InvalidParameterError, "posRefBegInd < 0");
473 }
474 if (static_cast<size_t>(posRefBegInd) >= posRefCat.size()) {
475 throw LSST_EXCEPT(pexExcept::InvalidParameterError, "posRefBegInd too big");
476 }
477 double const maxRotationRad = geom::degToRad(control.maxRotationDeg);
478
479 // Create an undistorted Wcs to project everything with
480 // We'll anchor it at the center of the area.
481 geom::Extent2D srcCenter(0, 0);
482 for (auto iter = sourceCat.begin(); iter != sourceCat.end(); ++iter) {
483 srcCenter += geom::Extent2D(iter->getCentroid());
484 }
485 srcCenter /= sourceCat.size();
486
487 sphgeom::Vector3d refCenter(0, 0, 0);
488 for (auto iter = posRefCat.begin(); iter != posRefCat.end(); ++iter) {
489 refCenter += iter->getCoord().getVector();
490 }
491 refCenter /= posRefCat.size();
492
493 auto tanWcs =
494 afw::geom::makeSkyWcs(geom::Point2D(srcCenter), geom::SpherePoint(refCenter), wcs.getCdMatrix());
495
496 ProxyVector posRefProxyCat = makeProxies(posRefCat, *tanWcs);
497 ProxyVector sourceProxyCat = makeProxies(sourceCat, wcs, *tanWcs);
498
499 // sourceSubCat contains at most the numBrightStars brightest sources, sorted by decreasing flux
500 ProxyVector sourceSubCat =
501 selectPoint(sourceProxyCat, sourceCat.getSchema().find<double>(control.sourceFluxField).key,
502 control.numBrightStars);
503
504 // posRefSubCat skips the initial posRefBegInd brightest reference objects and contains
505 // at most the next len(sourceSubCat) + 25 brightest reference objects, sorted by decreasing flux
506 ProxyVector posRefSubCat =
507 selectPoint(posRefProxyCat, posRefCat.getSchema().find<double>(control.refFluxField).key,
508 sourceSubCat.size() + 25, posRefBegInd);
509 if (verbose) {
510 std::cout << "Catalog sizes: " << sourceSubCat.size() << " " << posRefSubCat.size() << std::endl;
511 }
512
513 // Construct a list of pairs of position reference stars sorted by increasing separation
514 std::vector<ProxyPair> posRefPairList;
515 size_t const posRefCatSubSize = posRefSubCat.size();
516 for (size_t i = 0; i < posRefCatSubSize - 1; i++) {
517 for (size_t j = i + 1; j < posRefCatSubSize; j++) {
518 posRefPairList.push_back(ProxyPair(posRefSubCat[i], posRefSubCat[j]));
519 }
520 }
521 std::sort(posRefPairList.begin(), posRefPairList.end(), cmpPair);
522
523 // Construct a list of pairs of sources sorted by increasing separation
524 std::vector<ProxyPair> sourcePairList;
525 size_t const sourceSubCatSize = sourceSubCat.size();
526 for (size_t i = 0; i < sourceSubCatSize - 1; i++) {
527 for (size_t j = i + 1; j < sourceSubCatSize; j++) {
528 sourcePairList.push_back(ProxyPair(sourceSubCat[i], sourceSubCat[j]));
529 }
530 }
531 std::sort(sourcePairList.begin(), sourcePairList.end(), cmpPair);
532
536
537 size_t const fullShapeSize = control.numPointsForShape - 1; // Max size of shape array
538 for (size_t ii = 0; ii < sourcePairList.size(); ii++) {
539 ProxyPair p = sourcePairList[ii];
540
542 searchPair(posRefPairList, p, control.matchingAllowancePix, maxRotationRad);
543
544 // If candidate pairs are found
545 if (q.size() != 0) {
546 std::vector<ProxyPair> srcMatPair;
547 std::vector<ProxyPair> catMatPair;
548
549 // Go through candidate pairs
550 for (size_t l = 0; l < q.size(); l++) {
551 // sign matters, so don't use deltaAng; no need to wrap because
552 // the result is used with deltaAng later
553 double dpa = p.pa - q[l].pa;
554
555 srcMatPair.clear();
556 catMatPair.clear();
557
558 srcMatPair.push_back(p);
559 catMatPair.push_back(q[l]);
560
561 if (verbose) {
562 std::cout << "p dist: " << p.distance << " pa: " << p.pa << std::endl;
563 std::cout << "q dist: " << q[l].distance << " pa: " << q[l].pa << std::endl;
564 }
565
566 for (size_t k = 0; k < sourceSubCat.size(); k++) {
567 if (p.first == sourceSubCat[k] || p.second == sourceSubCat[k]) continue;
568
569 ProxyPair pp(p.first, sourceSubCat[k]);
570
571 std::vector<ProxyPair>::iterator r = searchPair3(
572 posRefPairList, pp, q[l], control.matchingAllowancePix, dpa, maxRotationRad);
573 if (r != posRefPairList.end()) {
574 srcMatPair.push_back(pp);
575 catMatPair.push_back(*r);
576 if (verbose) {
577 std::cout << " p dist: " << pp.distance << " pa: " << pp.pa << std::endl;
578 std::cout << " r dist: " << (*r).distance << " pa: " << (*r).pa << std::endl;
579 }
580 if (srcMatPair.size() == fullShapeSize) {
581 break;
582 }
583 }
584 }
585
586 bool goodMatch = false;
587 if (srcMatPair.size() == fullShapeSize) {
588 goodMatch = true;
589 for (size_t k = 1; k < catMatPair.size(); k++) {
590 if (catMatPair[0].first != catMatPair[k].first) {
591 goodMatch = false;
592 break;
593 }
594 }
595 }
596
597 if (goodMatch && srcMatPair.size() == fullShapeSize) {
598 ProxyVector srcMat;
599 ProxyVector catMat;
600
601 srcMat.push_back(srcMatPair[0].first);
602 catMat.push_back(catMatPair[0].first);
603 for (size_t k = 0; k < srcMatPair.size(); k++) {
604 srcMat.push_back(srcMatPair[k].second);
605 catMat.push_back(catMatPair[k].second);
606 }
607
608 boost::shared_array<double> coeff = polyfit(1, srcMat, catMat);
609
610 if (verbose) {
611 for (size_t k = 0; k < srcMat.size(); k++) {
612 std::cout << "circle(" << srcMat[k].getX() << "," << srcMat[k].getY()
613 << ",10) # color=green" << std::endl;
614 std::cout << "circle(" << catMat[k].getX() << "," << catMat[k].getY()
615 << ",10) # color=red" << std::endl;
616 std::cout << "line(" << srcMat[0].getX() << "," << srcMat[0].getY() << ","
617 << srcMat[k].getX() << "," << srcMat[k].getY()
618 << ") # line=0 0 color=green" << std::endl;
619 std::cout << "line(" << catMat[0].getX() << "," << catMat[0].getY() << ","
620 << catMat[k].getX() << "," << catMat[k].getY()
621 << ") # line=0 0 color=red" << std::endl;
622 }
623 }
624
625 double a = coeff[1];
626 double b = coeff[2];
627 double c = coeff[4];
628 double d = coeff[5];
629 geom::Angle const theta(std::acos((a * b + c * d) /
630 (std::sqrt(a * a + c * c) * std::sqrt(b * b + d * d))),
632 if (verbose) {
633 std::cout << "Linear fit from match:" << std::endl;
634 std::cout << coeff[0] << " " << coeff[1] << " " << coeff[2] << std::endl;
635 std::cout << coeff[3] << " " << coeff[4] << " " << coeff[5] << std::endl;
636 std::cout << "Determinant (max " << control.maxDeterminant << "): ";
637 std::cout << coeff[1] * coeff[5] - coeff[2] * coeff[4] - 1. << std::endl;
638 std::cout << "Angle between axes (deg; allowed 90 +/- ";
639 std::cout << control.allowedNonperpDeg << "): ";
640 std::cout << theta.asDegrees() << std::endl;
641 }
642 if (std::fabs(coeff[1] * coeff[5] - coeff[2] * coeff[4] - 1.) > control.maxDeterminant ||
643 std::fabs(theta.asDegrees() - 90) > control.allowedNonperpDeg ||
644 std::fabs(coeff[0]) > control.maxOffsetPix ||
645 std::fabs(coeff[3]) > control.maxOffsetPix) {
646 if (verbose) {
647 std::cout << "Bad; continuing" << std::endl;
648 }
649 continue;
650 } else {
651 double x0, y0, x1, y1;
652 int num = 0;
653 srcMat.clear();
654 catMat.clear();
655 for (size_t i = 0; i < sourceSubCat.size(); i++) {
656 x0 = sourceSubCat[i].getX();
657 y0 = sourceSubCat[i].getY();
658 transform(1, coeff, x0, y0, &x1, &y1);
659 auto refObjDist =
660 searchNearestPoint(posRefSubCat, x1, y1, control.matchingAllowancePix);
661 if (refObjDist.first != posRefSubCat.end()) {
662 num++;
663 srcMat.push_back(sourceSubCat[i]);
664 catMat.push_back(*refObjDist.first);
665 if (verbose) {
666 std::cout << "Match: " << x0 << "," << y0 << " --> " << x1 << "," << y1
667 << " <==> " << refObjDist.first->getX() << ","
668 << refObjDist.first->getY() << std::endl;
669 }
670 }
671 }
672 if (num <= control.numPointsForShape) {
673 // Can get matrix = 0,0,0,0; everything matches a single catalog object
674 if (verbose) {
675 std::cout << "Insufficient initial matches; continuing" << std::endl;
676 }
677 continue;
678 }
679 coeff = polyfit(1, srcMat, catMat);
680 if (verbose) {
681 std::cout << "Coefficients from initial matching:" << std::endl;
682 for (size_t i = 0; i < 6; ++i) {
683 std::cout << coeff[i] << " ";
684 }
686 }
687
688 matPair = FinalVerify(coeff, posRefProxyCat, sourceProxyCat,
689 control.matchingAllowancePix, verbose);
690 if (verbose) {
691 std::cout << "Number of matches: " << matPair.size() << " vs "
692 << control.minMatchedPairs << std::endl;
693 }
694 if (matPair.size() <= static_cast<std::size_t>(control.minMatchedPairs)) {
695 if (verbose) {
696 std::cout << "Insufficient final matches; continuing" << std::endl;
697 }
698 if (matPair.size() > matPairSave.size()) {
699 matPairSave = matPair;
700 }
701 continue;
702 } else {
703 if (verbose) {
704 std::cout << "Finish" << std::endl;
705 }
706 matPairCand.push_back(matPair);
707 if (matPairCand.size() == 3) {
708 goto END;
709 }
710 }
711 }
712 }
713 }
714 }
715 }
716
717END:
718 if (matPairCand.size() == 0) {
719 return matPairSave;
720 } else {
721 size_t nmatch = matPairCand[0].size();
722 afwTable::ReferenceMatchVector matPairRet = matPairCand[0];
723 for (size_t i = 1; i < matPairCand.size(); i++) {
724 if (matPairCand[i].size() > nmatch) {
725 nmatch = matPairCand[i].size();
726 matPairRet = matPairCand[i];
727 }
728 }
729 return matPairRet;
730 }
731}
table::Key< int > transform
table::Key< int > b
T acos(T... args)
A class representing an angle.
Definition Angle.h:128
Point in an unspecified spherical coordinate system.
Definition SpherePoint.h:57
T clear(T... args)
T endl(T... args)
T fabs(T... args)
Extent< double, 2 > Extent2D
Definition Extent.h:400
AngleUnit constexpr radians
constant with units of radians
Definition Angle.h:109
constexpr double degToRad(double x) noexcept
Definition Angle.h:52
ProxyVector makeProxies(afw::table::SourceCatalog const &sourceCat, afw::geom::SkyWcs const &distortedWcs, afw::geom::SkyWcs const &tanWcs)
T sort(T... args)
T sqrt(T... args)
table::Key< table::Array< double > > coeff
Definition PsfexPsf.cc:366

◆ PYBIND11_MODULE()

lsst::meas::astrom::PYBIND11_MODULE ( _measAstromLib ,
mod  )

Definition at line 25 of file _measAstromLib.cc.

25 {
26 WrapperCollection wrappers(mod, "lsst.meas.astrom");
27
28 wrappers.addInheritanceDependency("lsst.afw.math");
29
30 sip::wrapLeastSqFitter1d(wrappers);
31 sip::wrapCreateWcsWithSip(wrappers);
32 sip::wrapMatchSrcToCatalogue(wrappers);
33 sip::wrapLeastSqFitter2d(wrappers);
34
37 wrapSipTransform(wrappers);
38 wrapMatchOptimisticB(wrappers);
41 wrappers.finish();
42};
A helper class for subdividing pybind11 module across multiple translation units (i....
Definition python.h:242
void wrapSipTransform(WrapperCollection &wrappers)
void wrapMakeMatchStatistics(WrapperCollection &wrappers)
void wrapPolynomialTransform(WrapperCollection &wrappers)
void wrapScaledPolynomialTransformFitter(WrapperCollection &wrappers)
void wrapMatchOptimisticB(WrapperCollection &wrappers)
void wrapPessimisticPatternMatcherUtils(WrapperCollection &wrappers)

◆ rotateWcsPixelsBy90()

std::shared_ptr< afw::geom::SkyWcs > lsst::meas::astrom::rotateWcsPixelsBy90 ( afw::geom::SkyWcs const & wcs,
int nQuarter,
geom::Extent2I const & dimensions )

Return a new SkyWcs that represents a rotation of the image it corresponds to about the image's center.

Parameters
[in]wcsOriginal SkyWcs to be rotated.
[in]nQuarterNumber of 90 degree rotations (positive is counterclockwise).
[in]dimensionsWidth and height of the image.

Definition at line 179 of file SipTransform.cc.

180 {
181 geom::Extent2D offset;
182 switch (nQuarter % 4) {
183 case 0:
184 offset = geom::Extent2D(0, 0);
185 break;
186 case 1:
187 offset = geom::Extent2D(dimensions.getY() - 1, 0);
188 break;
189 case 2:
190 offset = geom::Extent2D(dimensions - geom::Extent2I(1, 1));
191 break;
192 case 3:
193 offset = geom::Extent2D(0, dimensions.getX() - 1);
194 break;
195 }
196 auto rot = geom::LinearTransform::makeRotation(nQuarter * 90.0 * geom::degrees);
197 return transformWcsPixels(wcs, geom::AffineTransform(rot, offset));
198}
AngleUnit constexpr degrees
constant with units of degrees
Definition Angle.h:110

◆ test_rotation()

RotationTestResult lsst::meas::astrom::test_rotation ( ndarray::Array< double, 1, 1 > const & src_center,
ndarray::Array< double, 1, 1 > const & ref_center,
ndarray::Array< double, 1, 1 > const & src_delta,
ndarray::Array< double, 1, 1 > const & ref_delta,
double cos_shift,
double max_cos_rot_sq )

Test if the rotation implied between the source pattern and reference pattern is within tolerance.

To test this we need to create the first part of our spherical rotation matrix which we also return for use later.

Parameters
[in]src_center3 vector defining the center of the candidate source pinwheel pattern.
[in]ref_center3 vector defining the center of the candidate reference pinwheel pattern.
[in]src_delta3 vector delta between the source pattern center and the end of the pinwheel spoke.
[in]ref_delta3 vector delta of the candidate matched reference pair
[in]cos_shiftCosine of the angle between the source and reference candidate centers.
[in]max_cos_rot_sqMaximum allowed rotation of the pinwheel pattern.
Returns
Result of the test rotation.

Definition at line 181 of file pessimisticPatternMatcherUtils.cc.

185 {
186 // Make sure the sine is a real number.
187 if (cos_shift > 1.0) {
188 cos_shift = 1;
189 } else if (cos_shift < -1.0) {
190 cos_shift = -1;
191 }
192 double sin_shift = sqrt(1 - cos_shift * cos_shift);
193
194 // If the sine of our shift is zero we only need to use the identity matrix for the shift. Else we
195 // construct the rotation matrix for shift.
196 auto ref_center_eigen = ndarray::asEigenMatrix(ref_center).head<3>();
197 Eigen::Matrix3d shift_matrix;
198 if (sin_shift > 0) {
199 Eigen::Vector3d rot_axis = ndarray::asEigenMatrix(src_center).head<3>().cross(ref_center_eigen);
200 rot_axis /= sin_shift;
201 shift_matrix = create_spherical_rotation_matrix(rot_axis, cos_shift, sin_shift);
202 } else {
203 shift_matrix = Eigen::Matrix3d::Identity();
204 }
205
206 // Now that we have our shift we apply it to the src delta vector and check the rotation.
207 Eigen::Vector3d rot_src_delta = shift_matrix * ndarray::asEigenMatrix(src_delta);
208 Eigen::Vector3d proj_src_delta = rot_src_delta - rot_src_delta.dot(ref_center_eigen) * ref_center_eigen;
209 Eigen::Vector3d ref_delta_eigen = ndarray::asEigenMatrix(ref_delta);
210
211 Eigen::Vector3d proj_ref_delta =
212 ref_delta_eigen - ref_delta_eigen.dot(ref_center_eigen) * ref_center_eigen;
213 double proj_src_delta_sq = proj_src_delta.dot(proj_ref_delta);
214 double cos_rot_sq = proj_src_delta_sq * proj_src_delta_sq /
215 (proj_src_delta.dot(proj_src_delta) * proj_ref_delta.dot(proj_ref_delta));
216 if (cos_rot_sq < max_cos_rot_sq) {
217 // Return failure if the rotation isn't within tolerance.
218 return RotationTestResult();
219 }
220 return RotationTestResult(cos_rot_sq, proj_ref_delta, shift_matrix);
221}

◆ transformWcsPixels()

std::shared_ptr< afw::geom::SkyWcs > lsst::meas::astrom::transformWcsPixels ( afw::geom::SkyWcs const & wcs,
geom::AffineTransform const & s )

Create a new SkyWcs whose pixel coordinate system has been transformed via an affine transform.

Parameters
[in]wcsOriginal SkyWcs object.
[in]sAffineTransform to apply to the pixel coordinate system.
Returns
a new Wcs that satisfies the following:
newWcs = transformWcsPixels(wcs, s);
assert(newWcs.skyToPixel(sky), s(wcs.skyToPixel(sky)));
assert(newWcs.pixelToSky(pixel), wcs.pixelToSky(s.inverted()(pixel)));
std::shared_ptr< afw::geom::SkyWcs > transformWcsPixels(afw::geom::SkyWcs const &wcs, geom::AffineTransform const &s)
Create a new SkyWcs whose pixel coordinate system has been transformed via an affine transform.
for all sky coordinates sky and pixel coordinates pixel.

Definition at line 173 of file SipTransform.cc.

174 {
175 auto affineTransform22 = afw::geom::makeTransform(s);
176 return afw::geom::makeModifiedWcs(*affineTransform22->inverted(), wcs, true);
177}

◆ wrapMakeMatchStatistics()

void lsst::meas::astrom::wrapMakeMatchStatistics ( WrapperCollection & wrappers)

Definition at line 48 of file makeMatchStatistics.cc.

48 {
49 auto &mod = wrappers.module;
50
51 declareMakeMatchStatistics<afw::table::ReferenceMatch>(mod);
52 declareMakeMatchStatistics<afw::table::SourceMatch>(mod);
53}

◆ wrapMatchOptimisticB()

void lsst::meas::astrom::wrapMatchOptimisticB ( WrapperCollection & wrappers)

Definition at line 95 of file matchOptimisticB.cc.

95 {
96 declareRecordProxy(wrappers);
97 declareProxyPair(wrappers);
98 declareMatchOptimisticBControl(wrappers);
99
100 wrappers.module.def("makeProxies",
102 afw::geom::SkyWcs const &)) &
103 makeProxies,
104 "sourceCat"_a, "distortedWcs"_a, "tanWcs"_a);
105 wrappers.module.def("makeProxies",
106 (ProxyVector(*)(afw::table::SimpleCatalog const &, afw::geom::SkyWcs const &)) & makeProxies,
107 "posRefCat"_a, "tanWcs"_a);
108
109 wrappers.module.def("matchOptimisticB", &matchOptimisticB, "posRefCat"_a, "sourceCat"_a, "control"_a, "wcs"_a,
110 "posRefBegInd"_a = 0, "verbose"_a = false);
111}
A 2-dimensional celestial WCS that transform pixels to ICRS RA/Dec, using the LSST standard for pixel...
Definition SkyWcs.h:117
Custom catalog class for record/table subclasses that are guaranteed to have an ID,...

◆ wrapPessimisticPatternMatcherUtils()

void lsst::meas::astrom::wrapPessimisticPatternMatcherUtils ( WrapperCollection & wrappers)

Definition at line 41 of file pessimisticPatternMatcherUtils.cc.

41 {
42 wrappers.wrapType(py::class_<PatternResult>(wrappers.module, "PatternResult"), [](auto &mod, auto &cls) {
43 cls.def_readonly("candidate_pairs", &PatternResult::candidate_pairs);
44 cls.def_readonly("shift_rot_matrix", &PatternResult::shift_rot_matrix);
45 cls.def_readonly("cos_shift", &PatternResult::cos_shift);
46 cls.def_readonly("sin_rot", &PatternResult::sin_rot);
47 cls.def_readonly("success", &PatternResult::success);
48 mod.def("construct_pattern_and_shift_rot_matrix", &construct_pattern_and_shift_rot_matrix,
49 "src_pattern_array"_a, "src_delta_array"_a, "src_dist_array"_a, "dist_array"_a, "id_array"_a,
50 "reference_array"_a, "n_match"_a, "max_cos_theta_shift"_a, "max_cos_rot_sq"_a, "max_dist_rad"_a);
51 });
52}

◆ wrapPolynomialTransform()

void lsst::meas::astrom::wrapPolynomialTransform ( WrapperCollection & wrappers)

Definition at line 106 of file polynomialTransform.cc.

106 {
107 declarePolynomialTransform(wrappers);
108 declareScaledPolynomialTransform(wrappers);
109
110 wrappers.module.def("compose",
111 (PolynomialTransform(*)(geom::AffineTransform const &, PolynomialTransform const &)) & compose,
112 "t1"_a, "t2"_a);
113 wrappers.module.def("compose",
114 (PolynomialTransform(*)(PolynomialTransform const &, geom::AffineTransform const &)) & compose,
115 "t1"_a, "t2"_a);
116}

◆ wrapScaledPolynomialTransformFitter()

void lsst::meas::astrom::wrapScaledPolynomialTransformFitter ( WrapperCollection & wrappers)

Definition at line 75 of file scaledPolynomialTransformFitter.cc.

75 {
76 declareOutlierRejectionControl(wrappers);
77 declareScaledPolynomialTransformFitter(wrappers);
78}

◆ wrapSipTransform()

void lsst::meas::astrom::wrapSipTransform ( WrapperCollection & wrappers)

Definition at line 111 of file sipTransform.cc.

111 {
112 declareSipTransformBase(wrappers);
113 declareSipForwardTransform(wrappers);
114 declareSipReverseTransform(wrappers);
115
116 wrappers.module.def("makeWcs", makeWcs, "sipForward"_a, "sipReverse"_a, "skyOrigin"_a);
117 wrappers.module.def("transformWcsPixels", transformWcsPixels, "wcs"_a, "s"_a);
118 wrappers.module.def("rotateWcsPixelsBy90", rotateWcsPixelsBy90, "wcs"_a, "nQuarter"_a, "dimensions"_a);
119}