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