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