lsst.jointcal  21.0.0-20-g3b2d1f0+e09d92ee06
AstrometryFit.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 /*
3  * This file is part of jointcal.
4  *
5  * Developed for the LSST Data Management System.
6  * This product includes software developed by the LSST Project
7  * (https://www.lsst.org).
8  * See the COPYRIGHT file at the top-level directory of this distribution
9  * for details of code ownership.
10  *
11  * This program is free software: you can redistribute it and/or modify
12  * it under the terms of the GNU General Public License as published by
13  * the Free Software Foundation, either version 3 of the License, or
14  * (at your option) any later version.
15  *
16  * This program is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19  * GNU General Public License for more details.
20  *
21  * You should have received a copy of the GNU General Public License
22  * along with this program. If not, see <https://www.gnu.org/licenses/>.
23  */
24 
25 #include <iostream>
26 #include <iomanip>
27 #include <algorithm>
28 #include <fstream>
29 
30 #include "Eigen/Sparse"
31 
32 #include "lsst/log/Log.h"
33 #include "lsst/pex/exceptions.h"
36 #include "lsst/jointcal/Chi2.h"
42 
43 namespace lsst {
44 namespace jointcal {
45 
47  std::shared_ptr<AstrometryModel> astrometryModel, double posError)
48  : FitterBase(associations),
49  _astrometryModel(astrometryModel),
50  _refractionCoefficient(0),
51  _nParRefrac(_associations->getNFilters()),
52  _posError(posError) {
53  _log = LOG_GET("jointcal.AstrometryFit");
54  _JDRef = 0;
55 
56  _posError = posError;
57 
58  _referenceColor = 0;
59  _sigCol = 0;
60  std::size_t count = 0;
61  for (auto const &i : _associations->fittedStarList) {
62  _referenceColor += i->color;
63  _sigCol += std::pow(i->color, 2);
64  count++;
65  }
66  if (count) {
67  _referenceColor /= double(count);
68  if (_sigCol > 0) _sigCol = sqrt(_sigCol / count - std::pow(_referenceColor, 2));
69  }
70  LOGLS_INFO(_log, "Reference Color: " << _referenceColor << " sig " << _sigCol);
71 }
72 
73 #define NPAR_PM 2
74 
75 /* ! this routine is used in 3 instances: when computing
76 the derivatives, when computing the Chi2, when filling a tuple.
77 */
78 Point AstrometryFit::transformFittedStar(FittedStar const &fittedStar, AstrometryTransform const &sky2TP,
79  Point const &refractionVector, double refractionCoeff,
80  double mjd) const {
81  Point fittedStarInTP = sky2TP.apply(fittedStar);
82  if (fittedStar.mightMove) {
83  fittedStarInTP.x += fittedStar.pmx * mjd;
84  fittedStarInTP.y += fittedStar.pmy * mjd;
85  }
86  // account for atmospheric refraction: does nothing if color
87  // have not been assigned
88  // the color definition shouldbe the same when computing derivatives
89  double color = fittedStar.color - _referenceColor;
90  fittedStarInTP.x += refractionVector.x * color * refractionCoeff;
91  fittedStarInTP.y += refractionVector.y * color * refractionCoeff;
92  return fittedStarInTP;
93 }
94 
98 static void tweakAstromMeasurementErrors(FatPoint &P, MeasuredStar const &Ms, double error) {
99  static bool called = false;
100  static double increment = 0;
101  if (!called) {
102  increment = std::pow(error, 2); // was in Preferences
103  called = true;
104  }
105  P.vx += increment;
106  P.vy += increment;
107 }
108 
109 // we could consider computing the chi2 here.
110 // (although it is not extremely useful)
111 void AstrometryFit::leastSquareDerivativesMeasurement(CcdImage const &ccdImage, TripletList &tripletList,
112  Eigen::VectorXd &fullGrad,
113  MeasuredStarList const *msList) const {
114  /**********************************************************************/
115  /* @note the math in this method and accumulateStatImage() must be kept consistent,
116  * in terms of +/- convention, definition of model, etc. */
117  /**********************************************************************/
118 
119  /* Setup */
120  /* this routine works in two different ways: either providing the
121  ccdImage, of providing the MeasuredStarList. In the latter case, the
122  ccdImage should match the one(s) in the list. */
123  if (msList) assert(&(msList->front()->getCcdImage()) == &ccdImage);
124 
125  // get the Mapping
126  const AstrometryMapping *mapping = _astrometryModel->getMapping(ccdImage);
127  // count parameters
128  std::size_t npar_mapping = (_fittingDistortions) ? mapping->getNpar() : 0;
129  std::size_t npar_pos = (_fittingPos) ? 2 : 0;
130  std::size_t npar_refrac = (_fittingRefrac) ? 1 : 0;
131  std::size_t npar_pm = (_fittingPM) ? NPAR_PM : 0;
132  std::size_t npar_tot = npar_mapping + npar_pos + npar_refrac + npar_pm;
133  // if (npar_tot == 0) this CcdImage does not contribute
134  // any constraint to the fit, so :
135  if (npar_tot == 0) return;
136  IndexVector indices(npar_tot, -1);
137  if (_fittingDistortions) mapping->getMappingIndices(indices);
138 
139  // proper motion stuff
140  double mjd = ccdImage.getMjd() - _JDRef;
141  // refraction stuff
142  Point refractionVector = ccdImage.getRefractionVector();
143  // transformation from sky to TP
144  auto sky2TP = _astrometryModel->getSkyToTangentPlane(ccdImage);
145  // reserve matrices once for all measurements
146  AstrometryTransformLinear dypdy;
147  // the shape of H (et al) is required this way in order to be able to
148  // separate derivatives along x and y as vectors.
149  Eigen::MatrixX2d H(npar_tot, 2), halpha(npar_tot, 2), HW(npar_tot, 2);
150  Eigen::Matrix2d transW(2, 2);
151  Eigen::Matrix2d alpha(2, 2);
152  Eigen::VectorXd grad(npar_tot);
153  // current position in the Jacobian
154  Eigen::Index kTriplets = tripletList.getNextFreeIndex();
155  const MeasuredStarList &catalog = (msList) ? *msList : ccdImage.getCatalogForFit();
156 
157  for (auto &i : catalog) {
158  const MeasuredStar &ms = *i;
159  if (!ms.isValid()) continue;
160  // tweak the measurement errors
161  FatPoint inPos = ms;
162  tweakAstromMeasurementErrors(inPos, ms, _posError);
163  H.setZero(); // we cannot be sure that all entries will be overwritten.
164  FatPoint outPos;
165  // should *not* fill H if whatToFit excludes mapping parameters.
166  if (_fittingDistortions)
167  mapping->computeTransformAndDerivatives(inPos, outPos, H);
168  else
169  mapping->transformPosAndErrors(inPos, outPos);
170 
171  std::size_t ipar = npar_mapping;
172  double det = outPos.vx * outPos.vy - std::pow(outPos.vxy, 2);
173  if (det <= 0 || outPos.vx <= 0 || outPos.vy <= 0) {
174  LOGLS_WARN(_log, "Inconsistent measurement errors: drop measurement at "
175  << Point(ms) << " in image " << ccdImage.getName());
176  continue;
177  }
178  transW(0, 0) = outPos.vy / det;
179  transW(1, 1) = outPos.vx / det;
180  transW(0, 1) = transW(1, 0) = -outPos.vxy / det;
181  // compute alpha, a triangular square root
182  // of transW (i.e. a Cholesky factor)
183  alpha(0, 0) = sqrt(transW(0, 0));
184  // checked that alpha*alphaT = transW
185  alpha(1, 0) = transW(0, 1) / alpha(0, 0);
186  // DB - I think that the next line is equivalent to : alpha(1,1) = 1./sqrt(outPos.vy)
187  // PA - seems correct !
188  alpha(1, 1) = 1. / sqrt(det * transW(0, 0));
189  alpha(0, 1) = 0;
190 
191  std::shared_ptr<FittedStar const> const fs = ms.getFittedStar();
192 
193  Point fittedStarInTP =
194  transformFittedStar(*fs, *sky2TP, refractionVector, _refractionCoefficient, mjd);
195 
196  // compute derivative of TP position w.r.t sky position ....
197  if (npar_pos > 0) // ... if actually fitting FittedStar position
198  {
199  sky2TP->computeDerivative(*fs, dypdy, 1e-3);
200  // sign checked
201  // TODO Still have to check with non trivial non-diagonal terms
202  H(npar_mapping, 0) = -dypdy.A11();
203  H(npar_mapping + 1, 0) = -dypdy.A12();
204  H(npar_mapping, 1) = -dypdy.A21();
205  H(npar_mapping + 1, 1) = -dypdy.A22();
206  indices[npar_mapping] = fs->getIndexInMatrix();
207  indices.at(npar_mapping + 1) = fs->getIndexInMatrix() + 1;
208  ipar += npar_pos;
209  }
210  /* only consider proper motions of objects allowed to move,
211  unless the fit is going to be degenerate */
212  if (_fittingPM && fs->mightMove) {
213  H(ipar, 0) = -mjd; // Sign unchecked but consistent with above
214  H(ipar + 1, 1) = -mjd;
215  indices[ipar] = fs->getIndexInMatrix() + 2;
216  indices[ipar + 1] = fs->getIndexInMatrix() + 3;
217  ipar += npar_pm;
218  }
219  if (_fittingRefrac) {
220  /* if the definition of color changes, it has to remain
221  consistent with transformFittedStar */
222  double color = fs->color - _referenceColor;
223  // sign checked
224  H(ipar, 0) = -refractionVector.x * color;
225  H(ipar, 1) = -refractionVector.y * color;
226  indices[ipar] = _refracPosInMatrix;
227  ipar += 1;
228  }
229 
230  // We can now compute the residual
231  Eigen::Vector2d res(fittedStarInTP.x - outPos.x, fittedStarInTP.y - outPos.y);
232 
233  // do not write grad = H*transW*res to avoid
234  // dynamic allocation of a temporary
235  halpha = H * alpha;
236  HW = H * transW;
237  grad = HW * res;
238  // now feed in triplets and fullGrad
239  for (std::size_t ipar = 0; ipar < npar_tot; ++ipar) {
240  for (std::size_t ic = 0; ic < 2; ++ic) {
241  double val = halpha(ipar, ic);
242  if (val == 0) continue;
243  tripletList.addTriplet(indices[ipar], kTriplets + ic, val);
244  }
245  fullGrad(indices[ipar]) += grad(ipar);
246  }
247  kTriplets += 2; // each measurement contributes 2 columns in the Jacobian
248  } // end loop on measurements
249  tripletList.setNextFreeIndex(kTriplets);
250 }
251 
252 void AstrometryFit::leastSquareDerivativesReference(FittedStarList const &fittedStarList,
253  TripletList &tripletList,
254  Eigen::VectorXd &fullGrad) const {
255  /**********************************************************************/
256  /* @note the math in this method and accumulateStatRefStars() must be kept consistent,
257  * in terms of +/- convention, definition of model, etc. */
258  /**********************************************************************/
259 
260  /* We compute here the derivatives of the terms involving fitted
261  stars and reference stars. They only provide contributions if we
262  are fitting positions: */
263  if (!_fittingPos) return;
264  /* the other case where the accumulation of derivatives stops
265  here is when there are no RefStars */
266  if (_associations->refStarList.size() == 0) return;
267  Eigen::Matrix2d W(2, 2);
268  Eigen::Matrix2d alpha(2, 2);
269  Eigen::Matrix2d H(2, 2), halpha(2, 2), HW(2, 2);
270  AstrometryTransformLinear der;
271  Eigen::Vector2d res, grad;
272  Eigen::Index indices[2 + NPAR_PM];
273  Eigen::Index kTriplets = tripletList.getNextFreeIndex();
274  /* We cannot use the spherical coordinates directly to evaluate
275  Euclidean distances, we have to use a projector on some plane in
276  order to express least squares. Not projecting could lead to a
277  disaster around the poles or across alpha=0. So we need a
278  projector. We construct a projector and will change its
279  projection point at every object */
280  TanRaDecToPixel proj(AstrometryTransformLinear(), Point(0., 0.));
281  for (auto const &i : fittedStarList) {
282  const FittedStar &fs = *i;
283  const RefStar *rs = fs.getRefStar();
284  if (rs == nullptr) continue;
285  proj.setTangentPoint(fs);
286  // fs projects to (0,0), no need to compute its transform.
287  FatPoint rsProj;
288  proj.transformPosAndErrors(*rs, rsProj);
289  // Compute the derivative of the projector to incorporate its effects on the errors.
290  proj.computeDerivative(fs, der, 1e-4);
291  // sign checked. TODO check that the off-diagonal terms are OK.
292  H(0, 0) = -der.A11();
293  H(1, 0) = -der.A12();
294  H(0, 1) = -der.A21();
295  H(1, 1) = -der.A22();
296  // TO DO : account for proper motions.
297  double det = rsProj.vx * rsProj.vy - std::pow(rsProj.vxy, 2);
298  if (rsProj.vx <= 0 || rsProj.vy <= 0 || det <= 0) {
299  LOGLS_WARN(_log, "RefStar error matrix not positive definite for: " << *rs);
300  continue;
301  }
302  W(0, 0) = rsProj.vy / det;
303  W(0, 1) = W(1, 0) = -rsProj.vxy / det;
304  W(1, 1) = rsProj.vx / det;
305  // compute alpha, a triangular square root
306  // of W (i.e. a Cholesky factor)
307  alpha(0, 0) = sqrt(W(0, 0));
308  // checked that alpha*alphaT = transW
309  alpha(1, 0) = W(0, 1) / alpha(0, 0);
310  alpha(1, 1) = 1. / sqrt(det * W(0, 0));
311  alpha(0, 1) = 0;
312  indices[0] = fs.getIndexInMatrix();
313  indices[1] = fs.getIndexInMatrix() + 1;
314  unsigned npar_tot = 2;
315  /* TODO: account here for proper motions in the reference
316  catalog. We can code the effect and set the value to 0. Most
317  (all?) catalogs do not even come with a reference epoch. Gaia
318  will change that. When refraction enters into the game, one should
319  pay attention to the orientation of the frame */
320 
321  /* The residual should be Proj(fs)-Proj(*rs) in order to be consistent
322  with the measurement terms. Since P(fs) = 0, we have: */
323  res[0] = -rsProj.x;
324  res[1] = -rsProj.y;
325  halpha = H * alpha;
326  // grad = H*W*res
327  HW = H * W;
328  grad = HW * res;
329  // now feed in triplets and fullGrad
330  for (std::size_t ipar = 0; ipar < npar_tot; ++ipar) {
331  for (unsigned ic = 0; ic < 2; ++ic) {
332  double val = halpha(ipar, ic);
333  if (val == 0) continue;
334  tripletList.addTriplet(indices[ipar], kTriplets + ic, val);
335  }
336  fullGrad(indices[ipar]) += grad(ipar);
337  }
338  kTriplets += 2; // each measurement contributes 2 columns in the Jacobian
339  }
340  tripletList.setNextFreeIndex(kTriplets);
341 }
342 
343 void AstrometryFit::accumulateStatImage(CcdImage const &ccdImage, Chi2Accumulator &accum) const {
344  /**********************************************************************/
347  /**********************************************************************/
348  /* Setup */
349  // 1 : get the Mapping's
350  const AstrometryMapping *mapping = _astrometryModel->getMapping(ccdImage);
351  // proper motion stuff
352  double mjd = ccdImage.getMjd() - _JDRef;
353  // refraction stuff
354  Point refractionVector = ccdImage.getRefractionVector();
355  // transformation from sky to TP
356  auto sky2TP = _astrometryModel->getSkyToTangentPlane(ccdImage);
357  // reserve matrix once for all measurements
358  Eigen::Matrix2Xd transW(2, 2);
359 
360  auto &catalog = ccdImage.getCatalogForFit();
361  for (auto const &ms : catalog) {
362  if (!ms->isValid()) continue;
363  // tweak the measurement errors
364  FatPoint inPos = *ms;
365  tweakAstromMeasurementErrors(inPos, *ms, _posError);
366 
367  FatPoint outPos;
368  // should *not* fill H if whatToFit excludes mapping parameters.
369  mapping->transformPosAndErrors(inPos, outPos);
370  double det = outPos.vx * outPos.vy - std::pow(outPos.vxy, 2);
371  if (det <= 0 || outPos.vx <= 0 || outPos.vy <= 0) {
372  LOGLS_WARN(_log, " Inconsistent measurement errors :drop measurement at "
373  << Point(*ms) << " in image " << ccdImage.getName());
374  continue;
375  }
376  transW(0, 0) = outPos.vy / det;
377  transW(1, 1) = outPos.vx / det;
378  transW(0, 1) = transW(1, 0) = -outPos.vxy / det;
379 
380  std::shared_ptr<FittedStar const> const fs = ms->getFittedStar();
381  Point fittedStarInTP =
382  transformFittedStar(*fs, *sky2TP, refractionVector, _refractionCoefficient, mjd);
383 
384  Eigen::Vector2d res(fittedStarInTP.x - outPos.x, fittedStarInTP.y - outPos.y);
385  double chi2Val = res.transpose() * transW * res;
386 
387  accum.addEntry(chi2Val, 2, ms);
388  } // end of loop on measurements
389 }
390 
391 void AstrometryFit::accumulateStatImageList(CcdImageList const &ccdImageList, Chi2Accumulator &accum) const {
392  for (auto const &ccdImage : ccdImageList) {
393  accumulateStatImage(*ccdImage, accum);
394  }
395 }
396 
397 void AstrometryFit::accumulateStatRefStars(Chi2Accumulator &accum) const {
398  /**********************************************************************/
401  /**********************************************************************/
402 
403  /* If you wonder why we project here, read comments in
404  AstrometryFit::leastSquareDerivativesReference(TripletList &TList, Eigen::VectorXd &Rhs) */
405  FittedStarList &fittedStarList = _associations->fittedStarList;
406  TanRaDecToPixel proj(AstrometryTransformLinear(), Point(0., 0.));
407  for (auto const &fs : fittedStarList) {
408  const RefStar *rs = fs->getRefStar();
409  if (rs == nullptr) continue;
410  proj.setTangentPoint(*fs);
411  // fs projects to (0,0), no need to compute its transform.
412  FatPoint rsProj;
413  proj.transformPosAndErrors(*rs, rsProj);
414  // TO DO : account for proper motions.
415  double rx = rsProj.x; // -fsProj.x (which is 0)
416  double ry = rsProj.y;
417  double det = rsProj.vx * rsProj.vy - std::pow(rsProj.vxy, 2);
418  double wxx = rsProj.vy / det;
419  double wyy = rsProj.vx / det;
420  double wxy = -rsProj.vxy / det;
421  double chi2 = wxx * std::pow(rx, 2) + 2 * wxy * rx * ry + wyy * std::pow(ry, 2);
422  accum.addEntry(chi2, 2, fs);
423  }
424 }
425 
427 
429 void AstrometryFit::getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, IndexVector &indices) const {
430  if (_fittingDistortions) {
431  const AstrometryMapping *mapping = _astrometryModel->getMapping(measuredStar.getCcdImage());
432  mapping->getMappingIndices(indices);
433  }
434  std::shared_ptr<FittedStar const> const fs = measuredStar.getFittedStar();
435  Eigen::Index fsIndex = fs->getIndexInMatrix();
436  if (_fittingPos) {
437  indices.push_back(fsIndex);
438  indices.push_back(fsIndex + 1);
439  }
440  // For securing the outlier removal, the next block is just useless
441  if (_fittingPM) {
442  for (std::size_t k = 0; k < NPAR_PM; ++k) indices.push_back(fsIndex + 2 + k);
443  }
444  /* Should not put the index of refaction stuff or we will not be
445  able to remove more than 1 star at a time. */
446 }
447 
449  _whatToFit = whatToFit;
450  LOGLS_INFO(_log, "assignIndices: Now fitting " << whatToFit);
451  _fittingDistortions = (_whatToFit.find("Distortions") != std::string::npos);
452  _fittingPos = (_whatToFit.find("Positions") != std::string::npos);
453  _fittingRefrac = (_whatToFit.find("Refrac") != std::string::npos);
454  if (_sigCol == 0 && _fittingRefrac) {
456  "Cannot fit refraction coefficients without a color lever arm. Ignoring refraction.");
457  _fittingRefrac = false;
458  }
459  _fittingPM = (_whatToFit.find("PM") != std::string::npos);
460  // When entering here, we assume that whatToFit has already been interpreted.
461 
462  _nModelParams = 0;
463  if (_fittingDistortions) _nModelParams = _astrometryModel->assignIndices(_whatToFit, 0);
464  std::size_t ipar = _nModelParams;
465 
466  if (_fittingPos) {
467  FittedStarList &fittedStarList = _associations->fittedStarList;
468  for (auto &fittedStar : fittedStarList) {
469  // the parameter layout here is used also
470  // - when filling the derivatives
471  // - when updating (offsetParams())
472  // - in GetMeasuredStarIndices
473  fittedStar->setIndexInMatrix(ipar);
474  ipar += 2;
475  if ((_fittingPM)&fittedStar->mightMove) ipar += NPAR_PM;
476  }
477  }
478  _nStarParams = ipar - _nModelParams;
479  if (_fittingRefrac) {
480  _refracPosInMatrix = ipar;
481  ipar += _nParRefrac;
482  }
483  _nTotal = ipar;
484  LOGLS_DEBUG(_log, "nParameters total: " << _nTotal << " model: " << _nModelParams
485  << " values: " << _nStarParams);
486 }
487 
488 void AstrometryFit::offsetParams(Eigen::VectorXd const &delta) {
489  if (delta.size() != _nTotal)
491  "AstrometryFit::offsetParams : the provided vector length is not compatible with "
492  "the current whatToFit setting");
493  if (_fittingDistortions) _astrometryModel->offsetParams(delta);
494 
495  if (_fittingPos) {
496  FittedStarList &fittedStarList = _associations->fittedStarList;
497  for (auto const &i : fittedStarList) {
498  FittedStar &fs = *i;
499  // the parameter layout here is used also
500  // - when filling the derivatives
501  // - when assigning indices (assignIndices())
502  Eigen::Index index = fs.getIndexInMatrix();
503  fs.x += delta(index);
504  fs.y += delta(index + 1);
505  if ((_fittingPM)&fs.mightMove) {
506  fs.pmx += delta(index + 2);
507  fs.pmy += delta(index + 3);
508  }
509  }
510  }
511  if (_fittingRefrac) {
512  _refractionCoefficient += delta(_refracPosInMatrix);
513  }
514 }
515 
517 #if (0)
518  const char *what2fit[] = {"Positions",
519  "Distortions",
520  "Refrac",
521  "Positions Distortions",
522  "Positions Refrac",
523  "Distortions Refrac",
524  "Positions Distortions Refrac"};
525 #endif
526  const char *what2fit[] = {"Positions", "Distortions", "Positions Distortions"};
527  // DEBUG
528  for (unsigned k = 0; k < sizeof(what2fit) / sizeof(what2fit[0]); ++k) {
529  assignIndices(what2fit[k]);
530  TripletList tripletList(10000);
531  Eigen::VectorXd grad(_nTotal);
532  grad.setZero();
533  leastSquareDerivatives(tripletList, grad);
534  SparseMatrixD jacobian(_nTotal, tripletList.getNextFreeIndex());
535  jacobian.setFromTriplets(tripletList.begin(), tripletList.end());
536  SparseMatrixD hessian = jacobian * jacobian.transpose();
537  LOGLS_DEBUG(_log, "npar : " << _nTotal << ' ' << _nModelParams);
538  }
539 }
540 
542  std::ofstream ofile(filename.c_str());
543  std::string separator = "\t";
544 
545  ofile << "#id" << separator << "xccd" << separator << "yccd " << separator;
546  ofile << "rx" << separator << "ry" << separator;
547  ofile << "xtp" << separator << "ytp" << separator;
548  ofile << "mag" << separator << "mjd" << separator;
549  ofile << "xErr" << separator << "yErr" << separator << "xyCov" << separator;
550  ofile << "xtpi" << separator << "ytpi" << separator;
551  ofile << "rxi" << separator << "ryi" << separator;
552  ofile << "color" << separator << "fsindex" << separator;
553  ofile << "ra" << separator << "dec" << separator;
554  ofile << "chi2" << separator << "nm" << separator;
555  ofile << "chip" << separator << "visit" << std::endl;
556 
557  ofile << "#id in source catalog" << separator << "coordinates in CCD (pixels)" << separator << separator;
558  ofile << "residual on TP (degrees)" << separator << separator;
559  ofile << "transformed coordinate in TP (degrees)" << separator << separator;
560  ofile << "rough magnitude" << separator << "Modified Julian Date of the measurement" << separator;
561  ofile << "transformed measurement uncertainty (degrees)" << separator << separator << separator;
562  ofile << "as-read position on TP (degrees)" << separator << separator;
563  ofile << "as-read residual on TP (degrees)" << separator << separator;
564  ofile << "currently unused" << separator << "unique index of the fittedStar" << separator;
565  ofile << "on sky position of fittedStar" << separator << separator;
566  ofile << "contribution to Chi2 (2D dofs)" << separator << "number of measurements of this fittedStar"
567  << separator;
568  ofile << "chip id" << separator << "visit id" << std::endl;
569 
570  const CcdImageList &ccdImageList = _associations->getCcdImageList();
571  for (auto const &ccdImage : ccdImageList) {
572  const MeasuredStarList &cat = ccdImage->getCatalogForFit();
573  const AstrometryMapping *mapping = _astrometryModel->getMapping(*ccdImage);
574  const auto readTransform = ccdImage->getReadWcs();
575  const Point &refractionVector = ccdImage->getRefractionVector();
576  double mjd = ccdImage->getMjd() - _JDRef;
577  for (auto const &ms : cat) {
578  if (!ms->isValid()) continue;
579  FatPoint tpPos;
580  FatPoint inPos = *ms;
581  tweakAstromMeasurementErrors(inPos, *ms, _posError);
582  mapping->transformPosAndErrors(inPos, tpPos);
583  auto sky2TP = _astrometryModel->getSkyToTangentPlane(*ccdImage);
584  const std::unique_ptr<AstrometryTransform> readPixToTangentPlane =
585  compose(*sky2TP, *readTransform);
586  FatPoint inputTpPos = readPixToTangentPlane->apply(inPos);
587  std::shared_ptr<FittedStar const> const fs = ms->getFittedStar();
588 
589  Point fittedStarInTP =
590  transformFittedStar(*fs, *sky2TP, refractionVector, _refractionCoefficient, mjd);
591  Point res = tpPos - fittedStarInTP;
592  Point inputRes = inputTpPos - fittedStarInTP;
593  double det = tpPos.vx * tpPos.vy - std::pow(tpPos.vxy, 2);
594  double wxx = tpPos.vy / det;
595  double wyy = tpPos.vx / det;
596  double wxy = -tpPos.vxy / det;
597  double chi2 = wxx * res.x * res.x + wyy * res.y * res.y + 2 * wxy * res.x * res.y;
598 
599  ofile << std::setprecision(9);
600  ofile << ms->getId() << separator << ms->x << separator << ms->y << separator;
601  ofile << res.x << separator << res.y << separator;
602  ofile << tpPos.x << separator << tpPos.y << separator;
603  ofile << fs->getMag() << separator << mjd << separator;
604  ofile << tpPos.vx << separator << tpPos.vy << separator << tpPos.vxy << separator;
605  ofile << inputTpPos.x << separator << inputTpPos.y << separator;
606  ofile << inputRes.x << separator << inputRes.y << separator;
607  ofile << fs->color << separator << fs->getIndexInMatrix() << separator;
608  ofile << fs->x << separator << fs->y << separator;
609  ofile << chi2 << separator << fs->getMeasurementCount() << separator;
610  ofile << ccdImage->getCcdId() << separator << ccdImage->getVisit() << std::endl;
611  } // loop on measurements in image
612  } // loop on images
613 }
614 
616  std::ofstream ofile(filename.c_str());
617  std::string separator = "\t";
618 
619  ofile << "#ra" << separator << "dec " << separator;
620  ofile << "rx" << separator << "ry" << separator;
621  ofile << "mag" << separator;
622  ofile << "xErr" << separator << "yErr" << separator << "xyCov" << separator;
623  ofile << "color" << separator << "fsindex" << separator;
624  ofile << "chi2" << separator << "nm" << std::endl;
625 
626  ofile << "#coordinates of fittedStar (degrees)" << separator << separator;
627  ofile << "residual on TP (degrees)" << separator << separator;
628  ofile << "magnitude" << separator;
629  ofile << "refStar transformed measurement uncertainty (degrees)" << separator << separator << separator;
630  ofile << "currently unused" << separator << "unique index of the fittedStar" << separator;
631  ofile << "refStar contribution to Chi2 (2D dofs)" << separator
632  << "number of measurements of this FittedStar" << std::endl;
633 
634  // The following loop is heavily inspired from AstrometryFit::computeChi2()
635  const FittedStarList &fittedStarList = _associations->fittedStarList;
637  for (auto const &i : fittedStarList) {
638  const FittedStar &fs = *i;
639  const RefStar *rs = fs.getRefStar();
640  if (rs == nullptr) continue;
641  proj.setTangentPoint(fs);
642  // fs projects to (0,0), no need to compute its transform.
643  FatPoint rsProj;
644  proj.transformPosAndErrors(*rs, rsProj);
645  double rx = rsProj.x; // -fsProj.x (which is 0)
646  double ry = rsProj.y;
647  double det = rsProj.vx * rsProj.vy - std::pow(rsProj.vxy, 2);
648  double wxx = rsProj.vy / det;
649  double wyy = rsProj.vx / det;
650  double wxy = -rsProj.vxy / det;
651  double chi2 = wxx * std::pow(rx, 2) + 2 * wxy * rx * ry + wyy * std::pow(ry, 2);
652 
653  ofile << std::setprecision(9);
654  ofile << fs.x << separator << fs.y << separator;
655  ofile << rx << separator << ry << separator;
656  ofile << fs.getMag() << separator;
657  ofile << rsProj.vx << separator << rsProj.vy << separator << rsProj.vxy << separator;
658  ofile << fs.color << separator << fs.getIndexInMatrix() << separator;
659  ofile << chi2 << separator << fs.getMeasurementCount() << std::endl;
660  } // loop on FittedStars
661 }
662 } // namespace jointcal
663 } // namespace lsst
#define NPAR_PM
ImageT val
Eigen::Matrix< double, Eigen::Dynamic, 2 > MatrixX2d
Definition: Eigenstuff.h:33
Eigen::SparseMatrix< double, 0, Eigen::Index > SparseMatrixD
Definition: Eigenstuff.h:35
#define LSST_EXCEPT(type,...)
#define LOGLS_WARN(logger, message)
#define LOG_GET(logger)
#define LOGLS_INFO(logger, message)
#define LOGLS_DEBUG(logger, message)
T begin(T... args)
T c_str(T... args)
void checkStuff()
DEBUGGING routine.
void saveChi2RefContributions(std::string const &filename) const override
Save a CSV file containing residuals of reference terms.
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the requested quantities.
AstrometryFit(std::shared_ptr< Associations > associations, std::shared_ptr< AstrometryModel > astrometryModel, double posError)
this is the only constructor
void assignIndices(std::string const &whatToFit) override
Set parameters to fit and assign indices in the big matrix.
void saveChi2MeasContributions(std::string const &filename) const override
Save a CSV file containing residuals of measurement terms.
virtual class needed in the abstraction of the distortion model
virtual void transformPosAndErrors(FatPoint const &where, FatPoint &outPoint) const =0
The same as above but without the parameter derivatives (used to evaluate chi^2)
implements the linear transformations (6 real coefficients).
double getMag() const
Definition: BaseStar.h:104
A Point with uncertainties.
Definition: FatPoint.h:34
The objects which have been measured several times.
Definition: FittedStar.h:61
const RefStar * getRefStar() const
Get the astrometric reference star associated with this star.
Definition: FittedStar.h:112
int getMeasurementCount() const
Definition: FittedStar.h:96
Eigen::Index getIndexInMatrix() const
Definition: FittedStar.h:106
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition: FittedStar.h:123
Base class for fitters.
Definition: FitterBase.h:53
void leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const
Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting.
Definition: FitterBase.cc:341
Eigen::Index _nModelParams
Definition: FitterBase.h:166
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:161
A list of MeasuredStar. They are usually filled in Associations::createCcdImage.
Definition: MeasuredStar.h:146
A point in a plane.
Definition: Point.h:37
double x
coordinate
Definition: Point.h:42
Objects used as position anchors, typically USNO stars.
Definition: RefStar.h:39
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane)
void transformPosAndErrors(const FatPoint &in, FatPoint &out) const
transform with analytical derivatives
void setTangentPoint(Point const &tangentPoint)
Resets the projection (or tangent) point.
Eigen::Index getNextFreeIndex() const
Definition: Tripletlist.h:47
T count(T... args)
T end(T... args)
T endl(T... args)
T find(T... args)
std::list< std::shared_ptr< CcdImage > > CcdImageList
Definition: CcdImage.h:46
std::unique_ptr< AstrometryTransform > compose(AstrometryTransform const &left, AstrometryTransform const &right)
Returns a pointer to a composition of transforms, representing left(right()).
Class for a simple mapping implementing a generic AstrometryTransform.
T pow(T... args)
T push_back(T... args)
T setprecision(T... args)
T sqrt(T... args)