lsst.meas.astrom  15.0-1-g1eca518
CreateWcsWithSip.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 
3 /*
4  * LSST Data Management System
5  * Copyright 2008, 2009, 2010 LSST Corporation.
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 #include <cmath>
25 
26 #include "Eigen/SVD"
27 #include "Eigen/Cholesky"
28 #include "Eigen/LU"
29 
32 #include "lsst/afw/geom/Angle.h"
34 #include "lsst/afw/geom/SkyWcs.h"
35 #include "lsst/log/Log.h"
37 
38 namespace {
39 LOG_LOGGER _log = LOG_GET("meas.astrom.sip");
40 }
41 
42 namespace lsst {
43 namespace meas {
44 namespace astrom {
45 namespace sip {
46 
47 namespace except = lsst::pex::exceptions;
48 namespace afwGeom = lsst::afw::geom;
49 namespace afwImg = lsst::afw::image;
50 namespace afwDet = lsst::afw::detection;
51 namespace afwMath = lsst::afw::math;
52 namespace afwTable = lsst::afw::table;
53 
54 namespace {
55 
56 double const MAX_DISTANCE_CRPIX_TO_BBOXCTR = 1000;
57 
58 /*
59  * Given an index and a SIP order, calculate p and q for the index'th term u^p v^q
60  * (Cf. Eqn 2 in http://fits.gsfc.nasa.gov/registry/sip/SIP_distortion_v1_0.pdf)
61  */
63 indexToPQ(int const index, int const order)
64 {
65  int p = 0, q = index;
66 
67  for (int decrement = order; q >= decrement && decrement > 0; --decrement) {
68  q -= decrement;
69  p++;
70  }
71 
72  return std::make_pair(p, q);
73 }
74 
75 Eigen::MatrixXd
76 calculateCMatrix(Eigen::VectorXd const& axis1, Eigen::VectorXd const& axis2, int const order)
77 {
78  int nTerms = 0.5*order*(order+1);
79 
80  int const n = axis1.size();
81  Eigen::MatrixXd C = Eigen::MatrixXd::Zero(n, nTerms);
82  for (int i = 0; i < n; ++i) {
83  for (int j = 0; j < nTerms; ++j) {
84  std::pair<int, int> pq = indexToPQ(j, order);
85  int p = pq.first, q = pq.second;
86 
87  assert(p + q < order);
88  C(i, j) = ::pow(axis1[i], p)*::pow(axis2[i], q);
89  }
90  }
91 
92  return C;
93 }
94 
101 Eigen::VectorXd
102 leastSquaresSolve(Eigen::VectorXd b, Eigen::MatrixXd A) {
103  assert(A.rows() == b.rows());
104  Eigen::VectorXd par = A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
105  return par;
106 }
107 
108 } // anonymous namespace
109 
111 template<class MatchT>
113  std::vector<MatchT> const & matches,
114  afw::geom::SkyWcs const& linearWcs,
115  int const order,
116  afwGeom::Box2I const& bbox,
117  int const ngrid
118 ):
119  _matches(matches),
120  _bbox(bbox),
121  _ngrid(ngrid),
122  _linearWcs(std::make_shared<afw::geom::SkyWcs>(linearWcs)),
123  _sipOrder(order+1),
124  _reverseSipOrder(order+2), //Higher order for reverse transform
125  _sipA(Eigen::MatrixXd::Zero(_sipOrder, _sipOrder)),
126  _sipB(Eigen::MatrixXd::Zero(_sipOrder, _sipOrder)),
127  _sipAp(Eigen::MatrixXd::Zero(_reverseSipOrder, _reverseSipOrder)),
128  _sipBp(Eigen::MatrixXd::Zero(_reverseSipOrder, _reverseSipOrder)),
129  _newWcs()
130 {
131  if (order < 2) {
132  throw LSST_EXCEPT(except::OutOfRangeError, "SIP must be at least 2nd order");
133  }
134  if (_sipOrder > 9) {
136  str(boost::format("SIP forward order %d exceeds the convention limit of 9") %
137  _sipOrder));
138  }
139  if (_reverseSipOrder > 9) {
141  str(boost::format("SIP reverse order %d exceeds the convention limit of 9") %
142  _reverseSipOrder));
143  }
144 
145  if (_matches.size() < std::size_t(_sipOrder)) {
146  throw LSST_EXCEPT(except::LengthError, "Number of matches less than requested sip order");
147  }
148 
149  if (_ngrid <= 0) {
150  _ngrid = 5*_sipOrder; // should be plenty
151  }
152 
153  /*
154  * We need a bounding box to define the region over which:
155  * The forward transformation should be valid
156  * We calculate the reverse transformartion
157  * If no BBox is provided, guess one from the input points (extrapolated a bit to allow for fact
158  * that a finite number of points won't reach to the edge of the image)
159  */
160  if (_bbox.isEmpty() && !_matches.empty() > 0) {
161  for (
162  typename std::vector<MatchT>::const_iterator ptr = _matches.begin();
163  ptr != _matches.end();
164  ++ptr
165  ) {
166  afwTable::SourceRecord const & src = *ptr->second;
167  _bbox.include(afwGeom::PointI(src.getX(), src.getY()));
168  }
169  float const borderFrac = 1/::sqrt(_matches.size()); // fractional border to add to exact BBox
170  afwGeom::Extent2I border(borderFrac*_bbox.getWidth(), borderFrac*_bbox.getHeight());
171 
172  _bbox.grow(border);
173  }
174 
175  // If crpix is too far from the center of the fit bbox, move it to the center to improve the fit
176  auto const initialCrpix = _linearWcs->getPixelOrigin();
177  auto const bboxCenter = afw::geom::Box2D(_bbox).getCenter();
178  if (std::hypot(initialCrpix[0] - bboxCenter[0], initialCrpix[1] - bboxCenter[1]) >
179  MAX_DISTANCE_CRPIX_TO_BBOXCTR) {
180  LOGL_DEBUG(_log,
181  "_linearWcs crpix = %d, %d farther than %f from bbox center; shifting to center %d, %d",
182  initialCrpix[0], initialCrpix[1], MAX_DISTANCE_CRPIX_TO_BBOXCTR, bboxCenter[0],
183  bboxCenter[1]);
184  _linearWcs = _linearWcs->getTanWcs(bboxCenter);
185  }
186 
187  // Calculate the forward part of the SIP distortion
188  _calculateForwardMatrices();
189 
190  // Build a new wcs incorporating the forward SIP matrix, it's all we know so far
191  auto const crval = _linearWcs->getSkyOrigin();
192  auto const crpix = _linearWcs->getPixelOrigin();
193  Eigen::MatrixXd cdMatrix = _linearWcs->getCdMatrix();
194  _newWcs = afw::geom::makeTanSipWcs(crpix, crval, cdMatrix, _sipA, _sipB);
195 
196  // Use _newWcs to calculate the forward transformation on a grid, and derive the back transformation
197  _calculateReverseMatrices();
198 
199  // Build a new wcs incorporating all SIP matrices
200  _newWcs = afw::geom::makeTanSipWcs(crpix, crval, cdMatrix, _sipA, _sipB, _sipAp, _sipBp);
201 }
202 
203 template<class MatchT>
204 void
206 {
207  // Assumes FITS (1-indexed) coordinates.
208  afwGeom::Point2D crpix = _linearWcs->getPixelOrigin();
209 
210  // Calculate u, v and intermediate world coordinates
211  int const nPoints = _matches.size();
212  Eigen::VectorXd u(nPoints), v(nPoints), iwc1(nPoints), iwc2(nPoints);
213 
214  int i = 0;
215  auto linearIwcToSky = getIntermediateWorldCoordsToSky(*_linearWcs);
216  for (
217  typename std::vector<MatchT>::const_iterator ptr = _matches.begin();
218  ptr != _matches.end();
219  ++ptr, ++i
220  ) {
221  afwTable::ReferenceMatch const & match = *ptr;
222 
223  // iwc: intermediate world coordinate positions of catalogue objects
224  auto c = match.first->getCoord();
225  afwGeom::Point2D p = linearIwcToSky->applyInverse(c);
226  iwc1[i] = p[0];
227  iwc2[i] = p[1];
228  // u and v are intermediate pixel coordinates of observed (distorted) positions
229  u[i] = match.second->getX() - crpix[0];
230  v[i] = match.second->getY() - crpix[1];
231  }
232  // Scale u and v down to [-1,,+1] in order to avoid too large numbers in the polynomials
233  double uMax = u.cwiseAbs().maxCoeff();
234  double vMax = v.cwiseAbs().maxCoeff();
235  double norm = std::max(uMax, vMax);
236  u = u/norm;
237  v = v/norm;
238 
239  // Forward transform
240  int ord = _sipOrder;
241  Eigen::MatrixXd forwardC = calculateCMatrix(u, v, ord);
242  Eigen::VectorXd mu = leastSquaresSolve(iwc1, forwardC);
243  Eigen::VectorXd nu = leastSquaresSolve(iwc2, forwardC);
244 
245  // Use mu and nu to refine CD
246 
247  // Given the implementation of indexToPQ(), the refined values
248  // of the elements of the CD matrices are in elements 1 and "_sipOrder" of mu and nu
249  // If the implementation of indexToPQ() changes, these assertions
250  // will catch that change.
251  assert ((indexToPQ(0, ord) == std::pair<int, int>(0, 0)));
252  assert ((indexToPQ(1, ord) == std::pair<int, int>(0, 1)));
253  assert ((indexToPQ(ord, ord) == std::pair<int, int>(1, 0)));
254 
255  // Scale back CD matrix
256  Eigen::Matrix2d CD;
257  CD(1,0) = nu[ord]/norm;
258  CD(1,1) = nu[1]/norm;
259  CD(0,0) = mu[ord]/norm;
260  CD(0,1) = mu[1]/norm;
261 
262  Eigen::Matrix2d CDinv = CD.inverse(); //Direct inverse OK for 2x2 matrix in Eigen
263 
264  // The zeroth elements correspond to a shift in crpix
265  crpix[0] -= mu[0]*CDinv(0,0) + nu[0]*CDinv(0,1);
266  crpix[1] -= mu[0]*CDinv(1,0) + nu[0]*CDinv(1,1);
267 
268  auto const crval = _linearWcs->getSkyOrigin();
269 
270  _linearWcs = afw::geom::makeSkyWcs(crpix, crval, CD);
271 
272  //Get Sip terms
273 
274  //The rest of the elements correspond to
275  //mu[i] == CD11*Apq + CD12*Bpq and
276  //nu[i] == CD21*Apq + CD22*Bpq and
277  //
278  //We solve for Apq and Bpq with the equation
279  // (Apq) = (CD11 CD12)-1 * (mu[i])
280  // (Bpq) (CD21 CD22) (nu[i])
281 
282  for(int i=1; i< mu.rows(); ++i) {
283  std::pair<int, int> pq = indexToPQ(i, ord);
284  int p = pq.first, q = pq.second;
285 
286  if(p + q > 1 && p + q < ord) {
287  Eigen::Vector2d munu(2,1);
288  munu(0) = mu(i);
289  munu(1) = nu(i);
290  Eigen::Vector2d AB = CDinv*munu;
291  // Scale back sip coefficients
292  _sipA(p,q) = AB[0]/::pow(norm,p+q);
293  _sipB(p,q) = AB[1]/::pow(norm,p+q);
294  }
295  }
296 }
297 
298 template<class MatchT>
300  int const ngrid2 = _ngrid*_ngrid;
301 
302  Eigen::VectorXd U(ngrid2), V(ngrid2);
303  Eigen::VectorXd delta1(ngrid2), delta2(ngrid2);
304 
305  int const x0 = _bbox.getMinX();
306  double const dx = _bbox.getWidth()/(double)(_ngrid - 1);
307  int const y0 = _bbox.getMinY();
308  double const dy = _bbox.getHeight()/(double)(_ngrid - 1);
309 
310  // wcs->getPixelOrigin() returns LSST-style (0-indexed) pixel coords.
311  afwGeom::Point2D crpix = _newWcs->getPixelOrigin();
312 
313  LOGL_DEBUG(_log, "_calcReverseMatrices: x0,y0 %i,%i, W,H %i,%i, ngrid %i, dx,dy %g,%g, CRPIX %g,%g",
314  x0, y0, _bbox.getWidth(), _bbox.getHeight(), _ngrid, dx, dy, crpix[0], crpix[1]);
315 
316  auto tanWcs = _newWcs->getTanWcs(_newWcs->getPixelOrigin());
317  auto applySipAB = afwGeom::makeWcsPairTransform(*_newWcs, *tanWcs);
318  int k = 0;
319  for (int i = 0; i < _ngrid; ++i) {
320  double const y = y0 + i*dy;
321  std::vector<afwGeom::Point2D> beforeSipABPoints;
322  beforeSipABPoints.reserve(_ngrid);
323  for (int j = 0; j < _ngrid; ++j) {
324  double const x = x0 + j*dx;
325  beforeSipABPoints.emplace_back(afwGeom::Point2D(x, y));
326  }
327  auto const afterSipABPoints = applySipAB->applyForward(beforeSipABPoints);
328  for (int j = 0; j < _ngrid; ++j, ++k) {
329  double const x = x0 + j*dx;
330  double u,v;
331  // u and v are intermediate pixel coordinates on a grid of positions
332  u = x - crpix[0];
333  v = y - crpix[1];
334 
335  // U and V are the result of applying the "forward" (A,B) SIP coefficients
336  afwGeom::Point2D const xy = afterSipABPoints[j];
337  U[k] = xy[0] - crpix[0];
338  V[k] = xy[1] - crpix[1];
339 
340  if ((i == 0 || i == (_ngrid-1) || i == (_ngrid/2)) &&
341  (j == 0 || j == (_ngrid-1) || j == (_ngrid/2))) {
342  LOGL_DEBUG(_log, " x,y (%.1f, %.1f), u,v (%.1f, %.1f), U,V (%.1f, %.1f)", x, y, u, v, U[k], V[k]);
343  }
344 
345  delta1[k] = u - U[k];
346  delta2[k] = v - V[k];
347  }
348  }
349 
350  // Scale down U and V in order to avoid too large numbers in the polynomials
351  double UMax = U.cwiseAbs().maxCoeff();
352  double VMax = V.cwiseAbs().maxCoeff();
353  double norm = (UMax > VMax) ? UMax : VMax;
354  U = U/norm;
355  V = V/norm;
356 
357  // Reverse transform
358  int const ord = _reverseSipOrder;
359  Eigen::MatrixXd reverseC = calculateCMatrix(U, V, ord);
360  Eigen::VectorXd tmpA = leastSquaresSolve(delta1, reverseC);
361  Eigen::VectorXd tmpB = leastSquaresSolve(delta2, reverseC);
362 
363  assert(tmpA.rows() == tmpB.rows());
364  for(int j=0; j< tmpA.rows(); ++j) {
365  std::pair<int, int> pq = indexToPQ(j, ord);
366  int p = pq.first, q = pq.second;
367  // Scale back sip coefficients
368  _sipAp(p, q) = tmpA[j]/::pow(norm,p+q);
369  _sipBp(p, q) = tmpB[j]/::pow(norm,p+q);
370  }
371 }
372 
373 template<class MatchT>
375  assert(_newWcs.get());
376  return makeMatchStatisticsInPixels(*_newWcs, _matches, afw::math::MEDIAN).getValue();
377 }
378 
379 template<class MatchT>
381  assert(_linearWcs.get());
382  return makeMatchStatisticsInPixels(*_linearWcs, _matches, afw::math::MEDIAN).getValue();
383 }
384 
385 template<class MatchT>
387  assert(_newWcs.get());
389  *_newWcs, _matches, afw::math::MEDIAN).getValue()*afw::geom::radians;
390 }
391 
392 template<class MatchT>
394  assert(_linearWcs.get());
396  *_linearWcs, _matches, afw::math::MEDIAN).getValue()*afw::geom::radians;
397 }
398 
399 
400 #define INSTANTIATE(MATCH) \
401  template class CreateWcsWithSip<MATCH>;
402 
405 
406 }}}}
407 
408 
double getScatterInPixels() const
Compute the median separation, in pixels, between items in this object&#39;s match list.
T empty(T... args)
AngleUnit constexpr radians
CreateWcsWithSip(std::vector< MatchT > const &matches, afw::geom::SkyWcs const &linearWcs, int const order, afw::geom::Box2I const &bbox=afw::geom::Box2I(), int const ngrid=0)
Construct a CreateWcsWithSip.
void include(Point2I const &point)
table::PointKey< double > crpix
T norm(const T &x)
STL namespace.
table::PointKey< double > crval
T end(T... args)
std::shared_ptr< Record2 > second
std::shared_ptr< RecordT > src
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.
Point2D const getCenter() const
T make_pair(T... args)
std::shared_ptr< Record1 > first
double x
std::shared_ptr< TransformPoint2ToSpherePoint > getIntermediateWorldCoordsToSky(SkyWcs const &wcs, bool simplify=true)
T max(T... args)
afw::geom::Angle getScatterOnSky() const
Compute the median on-sky separation between items in this object&#39;s match list.
double getValue(Property const prop=NOTHING) const
T hypot(T... args)
T size(T... args)
#define LSST_EXCEPT(type,...)
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.
void grow(int buffer)
T begin(T... args)
afw::geom::Angle getLinearScatterOnSky() const
Compute the median on-sky separation between items in this object&#39;s match list,.
std::shared_ptr< SkyWcs > makeSkyWcs(daf::base::PropertySet &metadata, bool strip=false)
Measure the distortions in an image plane and express them a SIP polynomials.
#define LOGL_DEBUG(logger, message...)
#define INSTANTIATE(MATCH)
#define LOG_GET(logger)
int y
std::shared_ptr< SkyWcs > makeTanSipWcs(Point2D const &crpix, SpherePoint const &crval, Eigen::Matrix2d const &cdMatrix, Eigen::MatrixXd const &sipA, Eigen::MatrixXd const &sipB)
T reserve(T... args)
double getLinearScatterInPixels() const
Compute the median radial separation between items in this object&#39;s match list.
T emplace_back(T... args)