lsst.meas.astrom  17.0.1-8-g7a3d54a+33
ScaledPolynomialTransformFitter.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 
3 /*
4  * LSST Data Management System
5  * Copyright 2016 LSST/AURA
6  *
7  * This product includes software developed by the
8  * LSST Project (http://www.lsst.org/).
9  *
10  * This program is free software: you can redistribute it and/or modify
11  * it under the terms of the GNU General Public License as published by
12  * the Free Software Foundation, either version 3 of the License, or
13  * (at your option) any later version.
14  *
15  * This program is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the LSST License Statement and
21  * the GNU General Public License along with this program. If not,
22  * see <http://www.lsstcorp.org/LegalNotices/>.
23  */
24 
25 #include <map>
26 
27 #include "Eigen/LU" // for determinant, even though it's a 2x2 that doesn't use actual LU implementation
28 
29 #include "boost/math/tools/minima.hpp"
30 
31 #include "lsst/geom/Box.h"
36 #include "lsst/afw/table/Match.h"
38 
39 // When developing this code, it was used to add an in-line check for the
40 // correctness of a particularly complex calculation that was difficult to
41 // factor out into something unit-testable (the larger context that code is in
42 // *is* unit tested, which should be guard against regressions). The test
43 // code is still in place, guarded by a check against this preprocessor macro;
44 // set it to 1 to re-enable that test during development, but be sure to
45 // set it back to zero before committing.
46 #define LSST_ScaledPolynomialTransformFitter_TEST_IN_PLACE 0
47 
48 namespace lsst {
49 namespace meas {
50 namespace astrom {
51 
52 // A singleton struct that manages the schema and keys for polynomial fitting catalogs.
54 public:
63  // We use uint16 instead of Flag since it's the only bool we have here, we
64  // may want NumPy views, and afw::table doesn't support [u]int8 fields.
66 
67  Keys(Keys const &) = delete;
68  Keys(Keys &&) = delete;
69  Keys &operator=(Keys const &) = delete;
70  Keys &operator=(Keys &&) = delete;
71 
72  static Keys const &forMatches() {
73  static Keys const it(0);
74  return it;
75  }
76 
77  static Keys const &forGrid() {
78  static Keys const it;
79  return it;
80  }
81 
82 private:
83  Keys(int)
84  : schema(),
85  refId(schema.addField<afw::table::RecordId>("ref_id", "ID of reference object in this match.")),
86  srcId(schema.addField<afw::table::RecordId>("src_id", "ID of source object in this match.")),
88  "source positions in pixel coordinates.", "pix")),
89  input(afw::table::Point2DKey::addFields(schema, "intermediate",
90  "reference positions in intermediate world coordinates",
91  "deg")),
93  schema, "initial", "reference positions transformed by initial WCS", "pix")),
95  schema, "model", "result of applying transform to reference positions", "pix")),
96  outputErr(
97  afw::table::CovarianceMatrixKey<float, 2>::addFields(schema, "src", {"x", "y"}, "pix")),
98  rejected(schema.addField<std::uint16_t>("rejected",
99  "True if the match should be rejected from the fit.")) {
100  schema.getCitizen().markPersistent();
101  }
102 
103  Keys()
104  : schema(),
106  schema, "output", "grid output positions in intermediate world coordinates", "deg")),
108  "grid input positions in pixel coordinates.", "pix")),
110  schema, "model", "result of applying transform to input positions", "deg")) {
111  schema.getCitizen().markPersistent();
112  }
113 };
114 
115 namespace {
116 
117 // Return the AffineTransforms that maps the given (x,y) coordinates to lie within (-1, 1)x(-1, 1)
120  for (auto const &record : data) {
121  bbox.include(geom::Point2D(record.get(key)));
122  };
123  return geom::AffineTransform(
124  geom::LinearTransform::makeScaling(0.5 * bbox.getWidth(), 0.5 * bbox.getHeight()))
125  .inverted() *
127 }
128 
129 } // namespace
130 
132  int maxOrder, afw::table::ReferenceMatchVector const &matches, afw::geom::SkyWcs const &initialWcs,
133  double intrinsicScatter) {
134  Keys const &keys = Keys::forMatches();
135  afw::table::BaseCatalog catalog(keys.schema);
136  catalog.reserve(matches.size());
137  float var2 = intrinsicScatter * intrinsicScatter;
138  auto initialIwcToSky = getIntermediateWorldCoordsToSky(initialWcs);
139  for (auto const &match : matches) {
140  auto record = catalog.addNew();
141  record->set(keys.refId, match.first->getId());
142  record->set(keys.srcId, match.second->getId());
143  record->set(keys.input, initialIwcToSky->applyInverse(match.first->getCoord()));
144  record->set(keys.initial, initialWcs.skyToPixel(match.first->getCoord()));
145  record->set(keys.output, match.second->getCentroid());
146  record->set(keys.outputErr, match.second->getCentroidErr() + var2 * Eigen::Matrix2f::Identity());
147  record->set(keys.rejected, false);
148  }
149  return ScaledPolynomialTransformFitter(catalog, keys, maxOrder, intrinsicScatter,
150  computeScaling(catalog, keys.input),
151  computeScaling(catalog, keys.output));
152 }
153 
155  int maxOrder, geom::Box2D const &bbox, int nGridX, int nGridY,
156  ScaledPolynomialTransform const &toInvert) {
157  Keys const &keys = Keys::forGrid();
158  afw::table::BaseCatalog catalog(keys.schema);
159  catalog.reserve(nGridX * nGridY);
160  geom::Extent2D dx(bbox.getWidth() / nGridX, 0.0);
161  geom::Extent2D dy(0.0, bbox.getHeight() / nGridY);
162  for (int iy = 0; iy < nGridY; ++iy) {
163  for (int ix = 0; ix < nGridX; ++ix) {
164  geom::Point2D point = bbox.getMin() + dx * ix + dy * iy;
165  auto record = catalog.addNew();
166  record->set(keys.output, point);
167  record->set(keys.input, toInvert(point));
168  }
169  }
170  return ScaledPolynomialTransformFitter(catalog, keys, maxOrder, 0.0, computeScaling(catalog, keys.input),
171  computeScaling(catalog, keys.output));
172 }
173 
174 ScaledPolynomialTransformFitter::ScaledPolynomialTransformFitter(afw::table::BaseCatalog const &data,
175  Keys const &keys, int maxOrder,
176  double intrinsicScatter,
177  geom::AffineTransform const &inputScaling,
178  geom::AffineTransform const &outputScaling)
179  : _keys(keys),
180  _intrinsicScatter(intrinsicScatter),
181  _data(data),
182  _outputScaling(outputScaling),
183  _transform(PolynomialTransform(maxOrder), inputScaling, outputScaling.inverted()),
184  _vandermonde(data.size(), detail::computePackedSize(maxOrder)) {
185  // Create a matrix that evaluates the max-order polynomials of all the (scaled) input positions;
186  // we'll extract subsets of this later when fitting to a subset of the matches and a lower order.
187  for (std::size_t i = 0; i < data.size(); ++i) {
188  geom::Point2D input = getInputScaling()(_data[i].get(_keys.input));
189  // x[k] == pow(x, k), y[k] == pow(y, k)
190  detail::computePowers(_transform._poly._u, input.getX());
191  detail::computePowers(_transform._poly._v, input.getY());
192  // We pack coefficients in the following order:
193  // (0,0), (0,1), (1,0), (0,2), (1,1), (2,0)
194  // Note that this lets us choose the just first N(N+1)/2 columns to
195  // evaluate an Nth order polynomial, even if N < maxOrder.
196  for (int n = 0, j = 0; n <= maxOrder; ++n) {
197  for (int p = 0, q = n; p <= n; ++p, --q, ++j) {
198  _vandermonde(i, j) = _transform._poly._u[p] * _transform._poly._v[q];
199  }
200  }
201  }
202 }
203 
205  int maxOrder = _transform.getPoly().getOrder();
206  if (order < 0) {
207  order = maxOrder;
208  }
209  if (order > maxOrder) {
210  throw LSST_EXCEPT(
212  (boost::format("Order (%d) exceeded maximum order for the fitter (%d)") % order % maxOrder)
213  .str());
214  }
215 
216  int const packedSize = detail::computePackedSize(order);
217  std::size_t nGood = 0;
218  if (_keys.rejected.isValid()) {
219  for (auto const &record : _data) {
220  if (!record.get(_keys.rejected)) {
221  ++nGood;
222  }
223  }
224  } else {
225  nGood = _data.size();
226  }
227  // One block of the block-diagonal (2x2) unweighted design matrix M;
228  // m[i,j] = u_i^{p(j)} v_i^{q(j)}. The two nonzero blocks are the same,
229  // because we're using the same polynomial basis for x and y.
230  Eigen::MatrixXd m = Eigen::MatrixXd::Zero(nGood, packedSize);
231  // vx, vy: (2x1) blocks of the unweighted data vector v
232  Eigen::VectorXd vx = Eigen::VectorXd::Zero(nGood);
233  Eigen::VectorXd vy = Eigen::VectorXd::Zero(nGood);
234  // sxx, syy, sxy: (2x2) blocks of the covariance matrix S, each of which is
235  // individually diagonal.
236  Eigen::ArrayXd sxx(nGood);
237  Eigen::ArrayXd syy(nGood);
238  Eigen::ArrayXd sxy(nGood);
239  Eigen::Matrix2d outS = _outputScaling.getLinear().getMatrix();
240  for (std::size_t i1 = 0, i2 = 0; i1 < _data.size(); ++i1) {
241  // check that the 'rejected' field (== 'not outlier-rejected') is both
242  // present in the schema and not rejected.
243  if (!_keys.rejected.isValid() || !_data[i1].get(_keys.rejected)) {
244  geom::Point2D output = _outputScaling(_data[i1].get(_keys.output));
245  vx[i2] = output.getX();
246  vy[i2] = output.getY();
247  m.row(i2) = _vandermonde.row(i1).head(packedSize);
248  if (_keys.outputErr.isValid()) {
249  Eigen::Matrix2d modelErr =
250  outS * _data[i1].get(_keys.outputErr).cast<double>() * outS.adjoint();
251  sxx[i2] = modelErr(0, 0);
252  sxy[i2] = modelErr(0, 1);
253  syy[i2] = modelErr(1, 1);
254  } else {
255  sxx[i2] = 1.0;
256  sxy[i2] = 0.0;
257  syy[i2] = 1.0;
258  }
259  ++i2;
260  }
261  }
262  // Do a blockwise inverse of S. Note that the result F is still symmetric
263  Eigen::ArrayXd fxx = 1.0 / (sxx - sxy.square() / syy);
264  Eigen::ArrayXd fyy = 1.0 / (syy - sxy.square() / sxx);
265  Eigen::ArrayXd fxy = -(sxy / sxx) * fyy;
266 #ifdef LSST_ScaledPolynomialTransformFitter_TEST_IN_PLACE
267  assert((sxx * fxx + sxy * fxy).isApproxToConstant(1.0));
268  assert((syy * fyy + sxy * fxy).isApproxToConstant(1.0));
269  assert((sxx * fxy).isApprox(-sxy * fyy));
270  assert((sxy * fxx).isApprox(-syy * fxy));
271 #endif
272  // Now that we've got all the block quantities, we'll form the full normal equations matrix.
273  // That's H = M^T F M:
274  Eigen::MatrixXd h(2 * packedSize, 2 * packedSize);
275  h.topLeftCorner(packedSize, packedSize) = m.adjoint() * fxx.matrix().asDiagonal() * m;
276  h.topRightCorner(packedSize, packedSize) = m.adjoint() * fxy.matrix().asDiagonal() * m;
277  h.bottomLeftCorner(packedSize, packedSize) = h.topRightCorner(packedSize, packedSize).adjoint();
278  h.bottomRightCorner(packedSize, packedSize) = m.adjoint() * fyy.matrix().asDiagonal() * m;
279  // And here's the corresponding RHS vector, g = M^T F v
280  Eigen::VectorXd g(2 * packedSize);
281  g.head(packedSize) = m.adjoint() * (fxx.matrix().asDiagonal() * vx + fxy.matrix().asDiagonal() * vy);
282  g.tail(packedSize) = m.adjoint() * (fxy.matrix().asDiagonal() * vx + fyy.matrix().asDiagonal() * vy);
283  // Solve the normal equations.
285  auto solution = lstsq.getSolution();
286  // Unpack the solution vector back into the polynomial coefficient matrices.
287  for (int n = 0, j = 0; n <= order; ++n) {
288  for (int p = 0, q = n; p <= n; ++p, --q, ++j) {
289  _transform._poly._xCoeffs(p, q) = solution[j];
290  _transform._poly._yCoeffs(p, q) = solution[j + packedSize];
291  }
292  }
293 }
294 
296  for (auto &record : _data) {
297  record.set(_keys.model, _transform(record.get(_keys.input)));
298  }
299 }
300 
302  if (!_keys.rejected.isValid()) {
304  "Cannot compute intrinsic scatter on fitter initialized with fromGrid.");
305  }
306  double newIntrinsicScatter = computeIntrinsicScatter();
307  float varDiff = newIntrinsicScatter * newIntrinsicScatter - _intrinsicScatter * _intrinsicScatter;
308  for (auto &record : _data) {
309  record.set(_keys.outputErr, record.get(_keys.outputErr) + varDiff * Eigen::Matrix2f::Identity());
310  }
311  _intrinsicScatter = newIntrinsicScatter;
312  return _intrinsicScatter;
313 }
314 
315 double ScaledPolynomialTransformFitter::computeIntrinsicScatter() const {
316  // We model the variance matrix of each match as the sum of the intrinsic variance (which we're
317  // trying to fit) and the per-source measurement uncertainty.
318  // We start by computing the variance directly, which yields the sum.
319  // At the same time, we find the maximum per-source variance. Since the per-source uncertainties
320  // are actually 2x2 matrices, we use the square of the major axis of that ellipse.
321  double directVariance = 0.0; // direct estimate of total scatter (includes measurement errors)
322  double maxMeasurementVariance = 0.0; // maximum of the per-match measurement uncertainties
323  double oldIntrinsicVariance = _intrinsicScatter * _intrinsicScatter;
324  std::size_t nGood = 0;
325  for (auto const &record : _data) {
326  if (!_keys.rejected.isValid() || !record.get(_keys.rejected)) {
327  auto delta = record.get(_keys.output) - record.get(_keys.model);
328  directVariance += 0.5 * delta.computeSquaredNorm();
329  double cxx = _keys.outputErr.getElement(record, 0, 0) - oldIntrinsicVariance;
330  double cyy = _keys.outputErr.getElement(record, 1, 1) - oldIntrinsicVariance;
331  double cxy = _keys.outputErr.getElement(record, 0, 1);
332  // square of semimajor axis of uncertainty error ellipse
333  double ca2 = 0.5 * (cxx + cyy + std::sqrt(cxx * cxx + cyy * cyy + 4 * cxy * cxy - 2 * cxx * cyy));
334  maxMeasurementVariance = std::max(maxMeasurementVariance, ca2);
335  ++nGood;
336  }
337  }
338  directVariance /= nGood;
339 
340  // Function that computes the -log likelihood of the current deltas with
341  // the variance modeled as described above.
342  auto logLikelihood = [this](double intrinsicVariance) {
343  // Uncertainties in the table right now include the old intrinsic scatter; need to
344  // subtract it off as we add the new one in.
345  double varDiff = intrinsicVariance - this->_intrinsicScatter * this->_intrinsicScatter;
346  double q = 0.0;
347  for (auto &record : this->_data) {
348  double x = record.get(this->_keys.output.getX()) - record.get(_keys.model.getX());
349  double y = record.get(this->_keys.output.getY()) - record.get(_keys.model.getY());
350  double cxx = this->_keys.outputErr.getElement(record, 0, 0) + varDiff;
351  double cyy = this->_keys.outputErr.getElement(record, 1, 1) + varDiff;
352  double cxy = this->_keys.outputErr.getElement(record, 0, 1);
353  double det = cxx * cyy - cxy * cxy;
354  q += (x * x * cyy - 2 * x * y * cxy + y * y * cxx) / det + std::log(det);
355  }
356  return q;
357  };
358 
359  // directVariance brackets the intrinsic variance from above, and this quantity
360  // brackets it from below:
361  double minIntrinsicVariance = std::max(0.0, directVariance - maxMeasurementVariance);
362 
363  // Minimize the negative log likelihood to find the best-fit intrinsic variance.
364  static constexpr int BITS_REQUIRED = 16; // solution good to ~1E-4
365  boost::uintmax_t maxIterations = 20;
366  auto result = boost::math::tools::brent_find_minima(logLikelihood, minIntrinsicVariance, directVariance,
367  BITS_REQUIRED, maxIterations);
368  return std::sqrt(result.first); // return RMS instead of variance
369 }
370 
372  OutlierRejectionControl const &ctrl) {
373  // If the 'rejected' field isn't present in the schema (because the fitter
374  // was constructed with fromGrid), we can't do outlier rejection.
375  if (!_keys.rejected.isValid()) {
377  "Cannot reject outliers on fitter initialized with fromGrid.");
378  }
379  if (static_cast<std::size_t>(ctrl.nClipMin) >= _data.size()) {
380  throw LSST_EXCEPT(
382  (boost::format("Not enough values (%d) to clip %d.") % _data.size() % ctrl.nClipMin).str());
383  }
385  for (auto &record : _data) {
386  Eigen::Matrix2d cov = record.get(_keys.outputErr).cast<double>();
387  Eigen::Vector2d d = (record.get(_keys.output) - record.get(_keys.model)).asEigen();
388  double r2 = d.dot(cov.inverse() * d);
389  rankings.insert(std::make_pair(r2, &record));
390  }
391  auto cutoff = rankings.upper_bound(ctrl.nSigma * ctrl.nSigma);
392  int nClip = 0, nGood = 0;
393  for (auto iter = rankings.begin(); iter != cutoff; ++iter) {
394  iter->second->set(_keys.rejected, false);
395  ++nGood;
396  }
397  for (auto iter = cutoff; iter != rankings.end(); ++iter) {
398  iter->second->set(_keys.rejected, true);
399  ++nClip;
400  }
401  assert(static_cast<std::size_t>(nGood + nClip) == _data.size());
402  while (nClip < ctrl.nClipMin) {
403  --cutoff;
404  cutoff->second->set(_keys.rejected, true);
405  ++nClip;
406  }
407  while (nClip > ctrl.nClipMax && cutoff != rankings.end()) {
408  cutoff->second->set(_keys.rejected, false);
409  ++cutoff;
410  --nClip;
411  }
413  if (cutoff != rankings.end()) {
414  result.first = std::sqrt(cutoff->first);
415  }
416  return result;
417 }
418 
419 } // namespace astrom
420 } // namespace meas
421 } // namespace lsst
daf::base::Citizen & getCitizen()
double getHeight() const noexcept
double updateIntrinsicScatter()
Infer the intrinsic scatter in the offset between the data points and the best-fit model...
double getWidth() const noexcept
int nClipMax
"Never clip more than this many matches." ;
lsst::geom::Point2D skyToPixel(lsst::geom::SpherePoint const &sky) const
Point2D const getCenter() const noexcept
static LeastSquares fromNormalEquations(ndarray::Array< T1, 2, C1 > const &fisher, ndarray::Array< T2, 1, C2 > const &rhs, Factorization factorization=NORMAL_EIGENSYSTEM)
int nClipMin
"Always clip at least this many matches." ;
T log(T... args)
table::Box2IKey bbox
void include(Point2D const &point) noexcept
static PointKey addFields(Schema &schema, std::string const &name, std::string const &doc, std::string const &unit)
std::pair< double, size_t > rejectOutliers(OutlierRejectionControl const &ctrl)
Mark outliers in the data catalog using sigma clipping.
geom::AffineTransform const & getInputScaling() const
Return the input scaling transform that maps input data points to [-1, 1].
Point2D const getMin() const noexcept
A fitter class for scaled polynomial transforms.
STL class.
void fit(int order=-1)
Perform a linear least-squares fit of the polynomial coefficients.
int computePackedSize(int order)
Compute this size of a packed 2-d polynomial coefficient array.
void updateModel()
Update the &#39;model&#39; field in the data catalog using the current best- fit transform.
static ScaledPolynomialTransformFitter fromMatches(int maxOrder, afw::table::ReferenceMatchVector const &matches, afw::geom::SkyWcs const &initialWcs, double intrinsicScatter)
Initialize a fit from intermediate world coordinates to pixels using source/reference matches...
char * data
def keys(self)
Key< U > key
T make_pair(T... args)
double x
std::shared_ptr< TransformPoint2ToSpherePoint > getIntermediateWorldCoordsToSky(SkyWcs const &wcs, bool simplify=true)
T max(T... args)
static LinearTransform makeScaling(double s) noexcept
std::unique_ptr< SchemaItem< U > > result
Eigen::Vector3d asEigen(sphgeom::Vector3d const &vector) noexcept
A 2-d coordinate transform represented by a lazy composition of an AffineTransform, a PolynomialTransform, and another AffineTransform.
void reserve(size_type n)
T size(T... args)
#define LSST_EXCEPT(type,...)
size_type size() const
T sqrt(T... args)
m
static ScaledPolynomialTransformFitter fromGrid(int maxOrder, geom::Box2D const &bbox, int nGridX, int nGridY, ScaledPolynomialTransform const &toInvert)
Initialize a fit that inverts an existing transform by evaluating and fitting to points on a grid...
Control object for outlier rejection in ScaledPolynomialTransformFitter.
int y
AffineTransform const inverted() const
Key< T > addField(Field< T > const &field, bool doReplace=false)
std::shared_ptr< BaseRecord > addNew()
void computePowers(Eigen::VectorXd &r, double x)
Fill an array with integer powers of x, so .
A 2-d coordinate transform represented by a pair of standard polynomials (one for each coordinate)...