lsst.meas.astrom  13.0-14-g9415442+26
 All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups Pages
matchOptimisticB.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <cmath>
3 #include <fstream>
4 #include <iostream>
5 #include <memory>
6 #include <utility>
7 #include <vector>
8 
9 #include "boost/shared_array.hpp"
10 #include "boost/multi_index_container.hpp"
11 #include "boost/multi_index/sequenced_index.hpp"
12 #include "boost/multi_index/ordered_index.hpp"
13 #include "boost/multi_index/global_fun.hpp"
14 
15 #include "gsl/gsl_linalg.h"
16 
17 #include "lsst/pex/exceptions.h"
18 #include "lsst/afw/image/Wcs.h"
19 #include "lsst/afw/image/DistortedTanWcs.h"
20 #include "lsst/afw/geom/Angle.h"
22 
23 namespace pexExcept = lsst::pex::exceptions;
24 namespace afwTable = lsst::afw::table;
25 
26 namespace {
27  using namespace lsst::meas::astrom;
28 
29  afwTable::RecordId getRefId(ProxyPair const &proxyPair) {
30  return proxyPair.first.record->getId();
31  }
32 
33  // a collection of ReferenceMatch indexed by insertion order and by reference object ID
34  struct refIdTag{}; // tag for retrieval by reference object ID
35  struct insertionOrderTag{}; // tag for retrieval by insertion order
36  typedef boost::multi_index_container<
37  ProxyPair,
38  boost::multi_index::indexed_by< // index by insertion order, so brightest sources first
39  boost::multi_index::sequenced<
40  boost::multi_index::tag<insertionOrderTag>
41  >,
42  boost::multi_index::ordered_unique< // index by reference object ID
43  boost::multi_index::tag<refIdTag>,
44  boost::multi_index::global_fun<
45  ProxyPair const &, afwTable::RecordId , &getRefId
46  > >
47  >
48  > MultiIndexedProxyPairList;
49 
50  // Algorithm is based on V.Tabur 2007, PASA, 24, 189-198
51  // "Fast Algorithms for Matching CCD Images to a Stellar Catalogue"
52 
59  inline double absDeltaAngle(double ang1, double ang2) {
60  return std::fmod(std::fabs(ang1 - ang2), M_PI*2);
61  }
62 
63  bool cmpPair(ProxyPair const &a, ProxyPair const &b) {
64  return a.distance > b.distance;
65  }
66 
72  struct CompareProxyFlux {
73 
74  bool operator()(RecordProxy const & a, RecordProxy const & b) const {
75  double aFlux = a.record->get(key);
76  double bFlux = b.record->get(key);
77  if (std::isnan(aFlux)) {
78  aFlux = 0.0;
79  }
80  if (std::isnan(bFlux)) {
81  bFlux = 0.0;
82  }
83  return aFlux > bFlux;
84  }
85 
86  afwTable::Key<double> key;
87  };
88 
99  ProxyVector selectPoint(
100  ProxyVector const &a,
101  afwTable::Key<double> const & key,
102  std::size_t num,
103  std::size_t startInd=0
104  ) {
105  if (startInd >= a.size()) {
106  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "startInd too big");
107  }
108  CompareProxyFlux cmp = {key};
109  ProxyVector b(a);
110  std::sort(b.begin(), b.end(), cmp);
111  std::size_t const endInd = std::min(startInd + num, b.size());
112  return ProxyVector(b.begin() + startInd, b.begin() + endInd);
113  }
114 
115  std::vector<ProxyPair> searchPair(
116  std::vector<ProxyPair> const &a,
117  ProxyPair const &p,
118  double e,
119  double e_dpa
120  ) {
121  std::vector<ProxyPair> v;
122 
123  for (size_t i = 0; i < a.size(); i++) {
124  double dd = std::fabs(a[i].distance - p.distance);
125  double dpa = absDeltaAngle(a[i].pa, p.pa);
126  if (dd < e && dpa < e_dpa) {
127  v.push_back(a[i]);
128  }
129  }
130 
131  return v;
132  }
133 
134  std::vector<ProxyPair>::iterator searchPair3(
135  std::vector<ProxyPair> &a,
136  ProxyPair const &p,
137  ProxyPair const &q,
138  double e,
139  double dpa,
140  double e_dpa = 0.02
141  ) {
142  std::vector<ProxyPair>::iterator idx = a.end();
143  double dd_min = 1.E+10;
144  //double dpa_min = e_dpa;
145 
146  for (auto i = a.begin(); i != a.end(); ++i) {
147  double dd = std::fabs(i->distance - p.distance);
148  #if 1
149  if (dd < e &&
150  absDeltaAngle(p.pa, i->pa - dpa) < e_dpa &&
151  dd < dd_min &&
152  (i->first == q.first)) {
153  dd_min = dd;
154  idx = i;
155  }
156  #else
157  if (dd < e &&
158  absDeltaAngle(p.pa, i->pa - dpa) < dpa_min) {
159  dpa_min = std::fabs(p.pa - i->pa - dpa);
160  idx = i;
161  }
162  #endif
163  }
164 
165  return idx;
166  }
167 
168  void transform(
169  int order,
170  boost::shared_array<double> const & coeff,
171  double x,
172  double y,
173  double *xn,
174  double *yn
175  ) {
176  int ncoeff = (order + 1) * (order + 2) / 2;
177  *xn = 0.0;
178  *yn = 0.0;
179  int n = 0;
180  for (int i = 0; i <= order; i++) {
181  for (int k = 0; k <= i; k++) {
182  int j = i - k;
183  *xn += coeff[n] * pow(x, j) * pow(y, k);
184  *yn += coeff[n+ncoeff] * pow(x, j) * pow(y, k);
185  n++;
186  }
187  }
188  }
189 
190  boost::shared_array<double> polyfit(
191  int order,
192  ProxyVector const &img,
193  ProxyVector const &posRefCat
194  ) {
195  int ncoeff = (order + 1) * (order + 2) / 2;
196  std::unique_ptr<int[]> xorder(new int[ncoeff]);
197  std::unique_ptr<int[]> yorder(new int[ncoeff]);
198 
199  int n = 0;
200  for (int i = 0; i <= order; i++) {
201  for (int k = 0; k <= i; k++) {
202  int j = i - k;
203  xorder[n] = j;
204  yorder[n] = k;
205  n++;
206  }
207  }
208 
209  std::unique_ptr<int[]> flag(new int[img.size()]);
210  for (size_t k = 0; k < img.size(); k++) {
211  flag[k] = 1;
212  }
213 
214  std::unique_ptr<double[]> a_data(new double[ncoeff*ncoeff]);
215  std::unique_ptr<double[]> b_data(new double[ncoeff]);
216  std::unique_ptr<double[]> c_data(new double[ncoeff]);
217 
218  boost::shared_array<double> coeff(new double[ncoeff*2]);
219 
220  for (int loop = 0; loop < 1; loop++) {
221  for (int i = 0; i < ncoeff; i++) {
222  for (int j = 0; j < ncoeff; j++) {
223  a_data[i*ncoeff+j] = 0.0;
224  for (size_t k = 0; k < img.size(); k++) {
225  if (flag[k] == 1) {
226  a_data[i*ncoeff+j] += pow(img[k].getX(), xorder[i]) *
227  pow(img[k].getY(), yorder[i]) *
228  pow(img[k].getX(), xorder[j]) *
229  pow(img[k].getY(), yorder[j]);
230  }
231  }
232  }
233  b_data[i] = c_data[i] = 0.0;
234  for (unsigned int k = 0; k < img.size(); k++) {
235  if (flag[k] == 1) {
236  b_data[i] += pow(img[k].getX(), xorder[i]) *
237  pow(img[k].getY(), yorder[i]) *
238  posRefCat[k].getX();
239  c_data[i] += pow(img[k].getX(), xorder[i]) *
240  pow(img[k].getY(), yorder[i]) *
241  posRefCat[k].getY();
242  }
243  }
244  }
245 
246  gsl_matrix_view a = gsl_matrix_view_array(a_data.get(), ncoeff, ncoeff);
247  gsl_vector_view b = gsl_vector_view_array(b_data.get(), ncoeff);
248  gsl_vector_view c = gsl_vector_view_array(c_data.get(), ncoeff);
249 
250  std::shared_ptr<gsl_vector> x(gsl_vector_alloc(ncoeff), gsl_vector_free);
251  std::shared_ptr<gsl_vector> y(gsl_vector_alloc(ncoeff), gsl_vector_free);
252 
253  int s;
254 
255  std::shared_ptr<gsl_permutation> p(gsl_permutation_alloc(ncoeff), gsl_permutation_free);
256 
257  gsl_linalg_LU_decomp(&a.matrix, p.get(), &s);
258  gsl_linalg_LU_solve(&a.matrix, p.get(), &b.vector, x.get());
259  gsl_linalg_LU_solve(&a.matrix, p.get(), &c.vector, y.get());
260 
261  for (int i = 0; i < ncoeff; i++) {
262  coeff[i] = x->data[i];
263  coeff[i+ncoeff] = y->data[i];
264  }
265 
266  double S, Sx, Sy, Sxx, Syy;
267  S = Sx = Sy = Sxx = Syy = 0.0;
268  for (size_t k = 0; k < img.size(); k++) {
269  if (flag[k] == 1) {
270  double x0 = img[k].getX();
271  double y0 = img[k].getY();
272  double x1, y1;
273  transform(order, coeff, x0, y0, &x1, &y1);
274  S += 1.;
275  Sx += (x1 - posRefCat[k].getX());
276  Sxx += (x1 - posRefCat[k].getX()) * (x1 - posRefCat[k].getX());
277  Sy += (y1 - posRefCat[k].getY());
278  Syy += (y1 - posRefCat[k].getY()) * (y1 - posRefCat[k].getY());
279  }
280  }
281  double x_sig = std::sqrt((Sxx - Sx * Sx / S) / S);
282  double y_sig = std::sqrt((Syy - Sy * Sy / S) / S);
283  //std::cout << x_sig << " " << y_sig << std::endl;
284 
285  for (size_t k = 0; k < img.size(); k++) {
286  double x0 = img[k].getX();
287  double y0 = img[k].getY();
288  double x1, y1;
289  transform(order, coeff, x0, y0, &x1, &y1);
290  if (std::fabs(x1-posRefCat[k].getX()) > 2. * x_sig ||
291  std::fabs(y1-posRefCat[k].getY()) > 2. * y_sig) {
292  flag[k] = 0;
293  }
294  }
295 
296  }
297 
298  return coeff;
299  }
300 
313  std::pair<ProxyVector::const_iterator, double> searchNearestPoint(
314  ProxyVector const &posRefCat,
315  double x,
316  double y,
317  double matchingAllowancePix
318  ) {
319  auto minDistSq = matchingAllowancePix*matchingAllowancePix;
320  auto foundPtr = posRefCat.end();
321  for (auto posRefPtr = posRefCat.begin(); posRefPtr != posRefCat.end(); ++posRefPtr) {
322  auto const dx = posRefPtr->getX() - x;
323  auto const dy = posRefPtr->getY() - y;
324  auto const distSq = dx*dx + dy*dy;
325  if (distSq < minDistSq) {
326  foundPtr = posRefPtr;
327  minDistSq = distSq;
328  }
329  }
330  return std::make_pair(foundPtr, std::sqrt(minDistSq));
331  }
332 
351  void addNearestMatch(
352  MultiIndexedProxyPairList &proxyPairList,
353  boost::shared_array<double> coeff,
354  ProxyVector const & posRefCat,
355  RecordProxy const &source,
356  double matchingAllowancePix
357  ) {
358  double x1, y1;
359  auto x0 = source.getX();
360  auto y0 = source.getY();
361  transform(1, coeff, x0, y0, &x1, &y1);
362  auto refObjDist = searchNearestPoint(posRefCat, x1, y1, matchingAllowancePix);
363  if (refObjDist.first == posRefCat.end()) {
364  // no reference object sufficiently close; do nothing
365  return;
366  }
367  auto existingMatch = proxyPairList.get<refIdTag>().find(refObjDist.first->record->getId());
368  if (existingMatch == proxyPairList.get<refIdTag>().end()) {
369  // reference object not used before; add new entry
370  auto proxyPair = ProxyPair(*refObjDist.first, source);
371  proxyPairList.get<refIdTag>().insert(proxyPair);
372  return;
373 #if 0
374  } else {
375  // This code keeps the closest match. It is disabled because it appears to be unnecessary
376  // and may be undesirable, since the image may have very faint sources which may give
377  // false matches. Note that the calling algorithm checks for matches starting with the
378  // brightest source and works down in source brightness.
379  if (existingMatch->distance <= refObjDist.second) {
380  // reference object used before and old source is closer; do nothing
381  return;
382  } else {
383  // reference object used before, but new source is closer; delete old entry and add new
384  proxyPairList.get<refIdTag>().erase(existingMatch);
385  auto proxyPair = ProxyPair(*refObjDist.first, source);
386  proxyPairList.get<refIdTag>().insert(proxyPair);
387  return;
388  }
389 #endif
390  }
391  }
392 
401  afwTable::ReferenceMatchVector FinalVerify(
402  boost::shared_array<double> coeff,
403  ProxyVector const & posRefCat,
404  ProxyVector const & sourceCat,
405  double matchingAllowancePix,
406  bool verbose
407  ) {
408  MultiIndexedProxyPairList proxyPairList;
409 
410  for (auto sourcePtr = sourceCat.begin(); sourcePtr != sourceCat.end(); ++sourcePtr) {
411  addNearestMatch(proxyPairList, coeff, posRefCat, *sourcePtr, matchingAllowancePix);
412  }
413  int order = 1;
414  if (proxyPairList.size() > 5) {
415  for (auto j = 0; j < 100; j++) {
416  auto prevNumMatches = proxyPairList.size();
417  auto srcMat = ProxyVector();
418  auto catMat = ProxyVector();
419  srcMat.reserve(proxyPairList.size());
420  catMat.reserve(proxyPairList.size());
421  for (auto matchPtr = proxyPairList.get<refIdTag>().begin();
422  matchPtr != proxyPairList.get<refIdTag>().end(); ++matchPtr) {
423  catMat.push_back(matchPtr->first);
424  srcMat.push_back(matchPtr->second);
425  }
426  coeff = polyfit(order, srcMat, catMat);
427  proxyPairList.clear();
428 
429  for (ProxyVector::const_iterator sourcePtr = sourceCat.begin();
430  sourcePtr != sourceCat.end(); ++sourcePtr) {
431  addNearestMatch(proxyPairList, coeff, posRefCat, *sourcePtr, matchingAllowancePix);
432  }
433  if (proxyPairList.size() == prevNumMatches) {
434  break;
435  }
436  }
437  }
438  auto matPair = afwTable::ReferenceMatchVector();
439  matPair.reserve(proxyPairList.size());
440  for (auto proxyPairIter = proxyPairList.get<insertionOrderTag>().begin();
441  proxyPairIter != proxyPairList.get<insertionOrderTag>().end(); ++proxyPairIter) {
442  matPair.push_back(afwTable::ReferenceMatch(
443  proxyPairIter->first.record,
444  std::static_pointer_cast<afwTable::SourceRecord>(proxyPairIter->second.record),
445  proxyPairIter->distance
446  ));
447  }
448 
449  return matPair;
450  }
451 
452 } // anonymous namespace
453 
454 namespace lsst {
455 namespace meas {
456 namespace astrom {
457 
458  void MatchOptimisticBControl::validate() const {
459  if (refFluxField.empty()) {
460  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "refFluxField must be specified");
461  }
462  if (sourceFluxField.empty()) {
463  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "sourceFluxField must be specified");
464  }
465  if (numBrightStars <= 0) {
466  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "numBrightStars must be positive");
467  }
468  if (minMatchedPairs < 0) {
469  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "minMatchedPairs must not be negative");
470  }
471  if (matchingAllowancePix <= 0) {
472  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "matchingAllowancePix must be positive");
473  }
474  if (maxOffsetPix <= 0) {
475  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "maxOffsetPix must be positive");
476  }
477  if (maxRotationDeg <= 0) {
478  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "maxRotationRad must be positive");
479  }
480  if (allowedNonperpDeg <= 0) {
481  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "allowedNonperpDeg must be positive");
482  }
483  if (numPointsForShape <= 0) {
484  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "numPointsForShape must be positive");
485  }
486  if (maxDeterminant <= 0) {
487  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "maxDeterminant must be positive");
488  }
489  }
490 
491  ProxyVector makeProxies(afwTable::SourceCatalog const & sourceCat,
492  afw::image::Wcs const& distortedWcs,
493  afw::image::Wcs const& tanWcs
494  )
495  {
496  ProxyVector r;
497  r.reserve(sourceCat.size());
498  for (afwTable::SourceCatalog::const_iterator sourcePtr = sourceCat.begin();
499  sourcePtr != sourceCat.end(); ++sourcePtr) {
500  r.push_back(RecordProxy(sourcePtr,
501  tanWcs.skyToPixel(*distortedWcs.pixelToSky(sourcePtr->getCentroid()))));
502  }
503  return r;
504  }
505 
506  ProxyVector makeProxies(afwTable::SimpleCatalog const & posRefCat,
507  afw::image::Wcs const& tanWcs
508  )
509  {
510  auto coordKey = afwTable::CoordKey(posRefCat.getSchema()["coord"]);
511  ProxyVector r;
512  r.reserve(posRefCat.size());
513  for (afwTable::SimpleCatalog::const_iterator posRefPtr = posRefCat.begin();
514  posRefPtr != posRefCat.end(); ++posRefPtr) {
515  r.push_back(RecordProxy(posRefPtr, tanWcs.skyToPixel(posRefPtr->get(coordKey))));
516  }
517  return r;
518  }
519 
520  afwTable::ReferenceMatchVector matchOptimisticB(
521  afwTable::SimpleCatalog const &posRefCat,
522  afwTable::SourceCatalog const &sourceCat,
523  MatchOptimisticBControl const &control,
524  afw::image::Wcs const& wcs,
525  int posRefBegInd,
526  bool verbose
527  ) {
528  control.validate();
529  if (posRefCat.empty()) {
530  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "no entries in posRefCat");
531  }
532  if (sourceCat.empty()) {
533  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "no entries in sourceCat");
534  }
535  if (posRefBegInd < 0) {
536  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "posRefBegInd < 0");
537  }
538  if (static_cast<size_t>(posRefBegInd) >= posRefCat.size()) {
539  throw LSST_EXCEPT(pexExcept::InvalidParameterError, "posRefBegInd too big");
540  }
541  double const maxRotationRad = afw::geom::degToRad(control.maxRotationDeg);
542 
543  CONST_PTR(afw::image::Wcs) tanWcs; // Undistorted Wcs, providing tangent plane
544  if (wcs.hasDistortion()) {
545  // Create an undistorted Wcs to project everything with
546  // We'll anchor it at the center of the area.
547  afw::geom::Extent2D srcCenter(0, 0);
548  for (auto iter = sourceCat.begin(); iter != sourceCat.end(); ++iter) {
549  srcCenter += afw::geom::Extent2D(iter->getCentroid());
550  }
551  srcCenter /= sourceCat.size();
552 
553  afw::geom::Extent3D refCenter(0, 0, 0);
554  for (auto iter = posRefCat.begin(); iter != posRefCat.end(); ++iter) {
555  refCenter += afw::geom::Extent3D(iter->getCoord().toIcrs().getVector());
556  }
557  refCenter /= posRefCat.size();
558 
559  tanWcs = std::make_shared<afw::image::Wcs>(
560  afw::coord::IcrsCoord(afw::geom::Point3D(refCenter)).getPosition(),
561  afw::geom::Point2D(srcCenter),
562  wcs.getCDMatrix()
563  );
564  }
565 
566  ProxyVector posRefProxyCat = makeProxies(posRefCat, tanWcs ? *tanWcs : wcs);
567  ProxyVector sourceProxyCat = makeProxies(sourceCat, wcs, tanWcs ? *tanWcs : wcs);
568 
569  // sourceSubCat contains at most the numBrightStars brightest sources, sorted by decreasing flux
570  ProxyVector sourceSubCat = selectPoint(
571  sourceProxyCat,
572  sourceCat.getSchema().find<double>(control.sourceFluxField).key,
573  control.numBrightStars);
574 
575  // posRefSubCat skips the initial posRefBegInd brightest reference objects and contains
576  // at most the next len(sourceSubCat) + 25 brightest reference objects, sorted by decreasing flux
577  ProxyVector posRefSubCat = selectPoint(
578  posRefProxyCat,
579  posRefCat.getSchema().find<double>(control.refFluxField).key,
580  sourceSubCat.size()+25,
581  posRefBegInd);
582  if (verbose) {
583  std::cout << "Catalog sizes: " << sourceSubCat.size() << " " << posRefSubCat.size() << std::endl;
584  }
585 
586  // Construct a list of pairs of position reference stars sorted by increasing separation
587  std::vector<ProxyPair> posRefPairList;
588  size_t const posRefCatSubSize = posRefSubCat.size();
589  for (size_t i = 0; i < posRefCatSubSize-1; i++) {
590  for (size_t j = i+1; j < posRefCatSubSize; j++) {
591  posRefPairList.push_back(ProxyPair(posRefSubCat[i], posRefSubCat[j]));
592  }
593  }
594  std::sort(posRefPairList.begin(), posRefPairList.end(), cmpPair);
595 
596  // Construct a list of pairs of sources sorted by increasing separation
597  std::vector<ProxyPair> sourcePairList;
598  size_t const sourceSubCatSize = sourceSubCat.size();
599  for (size_t i = 0; i < sourceSubCatSize-1; i++) {
600  for (size_t j = i+1; j < sourceSubCatSize; j++) {
601  sourcePairList.push_back(ProxyPair(sourceSubCat[i], sourceSubCat[j]));
602  }
603  }
604  std::sort(sourcePairList.begin(), sourcePairList.end(), cmpPair);
605 
606  afwTable::ReferenceMatchVector matPair;
607  afwTable::ReferenceMatchVector matPairSave;
608  std::vector<afwTable::ReferenceMatchVector> matPairCand;
609 
610  size_t const fullShapeSize = control.numPointsForShape - 1; // Max size of shape array
611  for (size_t ii = 0; ii < sourcePairList.size(); ii++) {
612  ProxyPair p = sourcePairList[ii];
613 
614  std::vector<ProxyPair> q = searchPair(posRefPairList, p,
615  control.matchingAllowancePix, maxRotationRad);
616 
617  // If candidate pairs are found
618  if (q.size() != 0) {
619 
620  std::vector<ProxyPair> srcMatPair;
621  std::vector<ProxyPair> catMatPair;
622 
623  // Go through candidate pairs
624  for (size_t l = 0; l < q.size(); l++) {
625  // sign matters, so don't use deltaAng; no need to wrap because
626  // the result is used with deltaAng later
627  double dpa = p.pa - q[l].pa;
628 
629  srcMatPair.clear();
630  catMatPair.clear();
631 
632  srcMatPair.push_back(p);
633  catMatPair.push_back(q[l]);
634 
635  if (verbose) {
636  std::cout << "p dist: " << p.distance << " pa: " << p.pa << std::endl;
637  std::cout << "q dist: " << q[l].distance << " pa: " << q[l].pa << std::endl;
638  }
639 
640  for (size_t k = 0; k < sourceSubCat.size(); k++) {
641  if (p.first == sourceSubCat[k] || p.second == sourceSubCat[k]) continue;
642 
643  ProxyPair pp(p.first, sourceSubCat[k]);
644 
645  std::vector<ProxyPair>::iterator r = searchPair3(posRefPairList, pp, q[l],
646  control.matchingAllowancePix, dpa, maxRotationRad);
647  if (r != posRefPairList.end()) {
648  srcMatPair.push_back(pp);
649  catMatPair.push_back(*r);
650  if (verbose) {
651  std::cout << " p dist: " << pp.distance << " pa: " << pp.pa << std::endl;
652  std::cout << " r dist: " << (*r).distance << " pa: " << (*r).pa << std::endl;
653  }
654  if (srcMatPair.size() == fullShapeSize) {
655  break;
656  }
657  }
658  }
659 
660  bool goodMatch = false;
661  if (srcMatPair.size() == fullShapeSize) {
662  goodMatch = true;
663  for (size_t k = 1; k < catMatPair.size(); k++) {
664  if (catMatPair[0].first != catMatPair[k].first) {
665  goodMatch = false;
666  break;
667  }
668  }
669  }
670 
671  if (goodMatch && srcMatPair.size() == fullShapeSize) {
672 
673  ProxyVector srcMat;
674  ProxyVector catMat;
675 
676  srcMat.push_back(srcMatPair[0].first);
677  catMat.push_back(catMatPair[0].first);
678  for (size_t k = 0; k < srcMatPair.size(); k++) {
679  srcMat.push_back(srcMatPair[k].second);
680  catMat.push_back(catMatPair[k].second);
681  }
682 
683  boost::shared_array<double> coeff = polyfit(1, srcMat, catMat);
684 
685  if (verbose) {
686  for (size_t k = 0; k < srcMat.size(); k++) {
687  std::cout << "circle(" << srcMat[k].getX() << ","
688  << srcMat[k].getY() << ",10) # color=green" << std::endl;
689  std::cout << "circle(" << catMat[k].getX() << ","
690  << catMat[k].getY() << ",10) # color=red" << std::endl;
691  std::cout << "line(" << srcMat[0].getX() << "," << srcMat[0].getY() << ","
692  << srcMat[k].getX() << "," << srcMat[k].getY()
693  << ") # line=0 0 color=green" << std::endl;
694  std::cout << "line(" << catMat[0].getX() << "," << catMat[0].getY() << ","
695  << catMat[k].getX() << "," << catMat[k].getY()
696  << ") # line=0 0 color=red" << std::endl;
697  }
698  }
699 
700  double a = coeff[1];
701  double b = coeff[2];
702  double c = coeff[4];
703  double d = coeff[5];
704  afw::geom::Angle const theta(std::acos((a*b+c*d)/
705  (std::sqrt(a*a+c*c)*std::sqrt(b*b+d*d))),
706  afw::geom::radians);
707  if (verbose) {
708  std::cout << "Linear fit from match:" << std::endl;
709  std::cout << coeff[0] << " " << coeff[1] << " " << coeff[2] << std::endl;
710  std::cout << coeff[3] << " " << coeff[4] << " " << coeff[5] << std::endl;
711  std::cout << "Determinant (max " << control.maxDeterminant << "): ";
712  std::cout << coeff[1] * coeff[5] - coeff[2] * coeff[4] - 1. << std::endl;
713  std::cout << "Angle between axes (deg; allowed 90 +/- ";
714  std::cout << control.allowedNonperpDeg << "): ";
715  std::cout << theta.asDegrees() << std::endl;
716  }
717  if (std::fabs(coeff[1] * coeff[5] - coeff[2] * coeff[4] - 1.)
718  > control.maxDeterminant ||
719  std::fabs(theta.asDegrees() - 90) > control.allowedNonperpDeg ||
720  std::fabs(coeff[0]) > control.maxOffsetPix ||
721  std::fabs(coeff[3]) > control.maxOffsetPix) {
722  if (verbose) {
723  std::cout << "Bad; continuing" << std::endl;
724  }
725  continue;
726  } else {
727  double x0, y0, x1, y1;
728  int num = 0;
729  srcMat.clear();
730  catMat.clear();
731  for (size_t i = 0; i < sourceSubCat.size(); i++) {
732  x0 = sourceSubCat[i].getX();
733  y0 = sourceSubCat[i].getY();
734  transform(1, coeff, x0, y0, &x1, &y1);
735  auto refObjDist = searchNearestPoint(posRefSubCat, x1, y1,
736  control.matchingAllowancePix);
737  if (refObjDist.first != posRefSubCat.end()) {
738  num++;
739  srcMat.push_back(sourceSubCat[i]);
740  catMat.push_back(*refObjDist.first);
741  if (verbose) {
742  std::cout << "Match: " << x0 << "," << y0 << " --> "
743  << x1 << "," << y1 << " <==> "
744  << refObjDist.first->getX() << "," << refObjDist.first->getY()
745  << std::endl;
746  }
747  }
748  }
749  if (num <= control.numPointsForShape) {
750  // Can get matrix = 0,0,0,0; everything matches a single catalog object
751  if (verbose) {
752  std::cout << "Insufficient initial matches; continuing" << std::endl;
753  }
754  continue;
755  }
756  coeff = polyfit(1, srcMat, catMat);
757  if (verbose) {
758  std::cout << "Coefficients from initial matching:" << std::endl;
759  for (size_t i = 0; i < 6; ++i) {
760  std::cout << coeff[i] << " ";
761  }
762  std::cout << std::endl;
763  }
764 
765  matPair = FinalVerify(coeff, posRefProxyCat, sourceProxyCat,
766  control.matchingAllowancePix, verbose);
767  if (verbose) {
768  std::cout << "Number of matches: " << matPair.size() << " vs " <<
769  control.minMatchedPairs << std::endl;
770  }
771  if (matPair.size() <= static_cast<std::size_t>(control.minMatchedPairs)) {
772  if (verbose) {
773  std::cout << "Insufficient final matches; continuing" << std::endl;
774  }
775  if (matPair.size() > matPairSave.size()) {
776  matPairSave = matPair;
777  }
778  continue;
779  } else {
780  if (verbose) {
781  std::cout << "Finish" << std::endl;
782  }
783  matPairCand.push_back(matPair);
784  if (matPairCand.size() == 3) {
785  goto END;
786  }
787  }
788  }
789  }
790  }
791  }
792  }
793 
794  END:
795  if (matPairCand.size() == 0) {
796  return matPairSave;
797  } else {
798  size_t nmatch = matPairCand[0].size();
799  afwTable::ReferenceMatchVector matPairRet = matPairCand[0];
800  for (size_t i = 1; i < matPairCand.size(); i++) {
801  if (matPairCand[i].size() > nmatch) {
802  nmatch = matPairCand[i].size();
803  matPairRet = matPairCand[i];
804  }
805  }
806  return matPairRet;
807  }
808  }
809 
810 }}} // namespace lsst::afw::astrom
int numPointsForShape
&quot;number of points in a matching shape&quot; ;
double matchingAllowancePix
&quot;maximum allowed distance between reference objects and sources (pixels)&quot; ;
std::vector< RecordProxy > ProxyVector
lsst::afw::table::ReferenceMatchVector matchOptimisticB(lsst::afw::table::SimpleCatalog const &posRefCat, lsst::afw::table::SourceCatalog const &sourceCat, MatchOptimisticBControl const &control, afw::image::Wcs const &wcs, int posRefBegInd=0, bool verbose=false)
Match sources to stars in a position reference catalog using optimistic pattern matching B...
double maxOffsetPix
&quot;maximum allowed frame translation (pixels)&quot; ;
int numBrightStars
&quot;maximum number of bright reference stars to use&quot; ;
std::string refFluxField
&quot;name of flux field in reference catalog&quot; ;
A wrapper around a SimpleRecord or SourceRecord that allows us to record a pixel position in a way th...
double maxRotationDeg
&quot;maximum allowed frame rotation (deg)&quot; ;
ProxyVector makeProxies(lsst::afw::table::SourceCatalog const &sourceCat, afw::image::Wcs const &distortedWcs, afw::image::Wcs const &tanWcs)
std::string sourceFluxField
&quot;name of flux field in source catalog&quot; ;
double allowedNonperpDeg
&quot;allowed non-perpendicularity of x and y axes (deg)&quot; ;
int minMatchedPairs
&quot;minimum number of matches&quot; ;
boost::shared_ptr< lsst::afw::table::SimpleRecord > record