lsst.jointcal  master-gc935ebf72c
 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 
5 #include "lsst/log/Log.h"
9 
10 #include "lsst/jointcal/Gtransfo.h"
11 #include "Eigen/Sparse"
12 #include "Eigen/CholmodSupport" // to switch to cholmod
13 #include "lsst/pex/exceptions.h"
14 #include <fstream>
16 
17 typedef Eigen::SparseMatrix<double> SpMat;
18 
19 #if 0
20 template<typename _MatrixType, int _UpLo = Eigen::Lower>
21 class CholmodSupernodalLLT2 : public Eigen::CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT2<_MatrixType, _UpLo> >
22 {
23  typedef Eigen::CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT2> Base;
24  using Base::m_cholmod;
25 
26 public:
27 
28  typedef _MatrixType MatrixType;
29  typedef typename MatrixType::Index Index;
30 
31  CholmodSupernodalLLT2() : Base() { init(); }
33  CholmodSupernodalLLT2(const MatrixType& matrix) : Base()
34  {
35  init();
36  this->compute(matrix);
37  }
38 
40  int update(const SpMat &H, const bool UpOrDown)
41  {
42  // check size
43  const Index size = Base::m_cholmodFactor->n;
44  EIGEN_UNUSED_VARIABLE(size);
45  eigen_assert(size == H.rows());
46  cholmod_sparse C_cs = viewAsCholmod(H);
47  /* We have to apply the magic permutation to the update matrix,
48  read page 117 of Cholmod UserGuide.pdf */
49  cholmod_sparse *C_cs_perm = cholmod_submatrix(C_cs,
50  Base::m_cholmodFactor->Perm,
51  Base::m_cholmodFactor->n,
52  nullptr, -1, nullptr, true, true,
53  &this->cholmod());
54 
55  int ret = cholmod_updown(UpOrDown, &C_cs_perm, Base::m_cholmodFactor, &this->cholmod());
56  cholmod_free_sparse(C_cs_perm, &this->cholmod());
57  return ret;
58  }
59 
60  ~CholmodSupernodalLLT2() {}
61 protected:
62  void init()
63  {
64  m_cholmod.final_asis = 1;
65  m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
66  // In CholmodBase::CholmodBase(), the following statement is missing in
67  // SuiteSparse 3.2.0.8. Fixed in 3.2.7
68  Base::m_shiftOffset[0] = Base::m_shiftOffset[1] = RealScalar(0.0);
69  }
70 };
71 #endif
72 
74 
78 template <typename _MatrixType, int _UpLo = Eigen::Lower>
80  : public Eigen::CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT2<_MatrixType, _UpLo> > {
81  typedef Eigen::CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT2> Base;
82  using Base::m_cholmod;
83 
84 public:
85  typedef _MatrixType MatrixType;
86  typedef typename MatrixType::Index Index;
87  typedef typename MatrixType::RealScalar RealScalar;
88 
89  CholmodSimplicialLDLT2() : Base() { init(); }
90 
91  CholmodSimplicialLDLT2(const MatrixType &matrix) : Base() {
92  init();
93  this->compute(matrix);
94  }
95 
96  // this routine is the one we added
97  int update(const SpMat &H, const bool UpOrDown) {
98  // check size
99  const Index size = Base::m_cholmodFactor->n;
100  EIGEN_UNUSED_VARIABLE(size);
101  eigen_assert(size == H.rows());
102 
103  cholmod_sparse C_cs = viewAsCholmod(H);
104  /* We have to apply the magic permutation to the update matrix,
105  read page 117 of Cholmod UserGuide.pdf */
106  cholmod_sparse *C_cs_perm =
107  cholmod_submatrix(&C_cs, (int *)Base::m_cholmodFactor->Perm, Base::m_cholmodFactor->n,
108  nullptr, -1, true, true, &this->cholmod());
109  assert(C_cs_perm);
110  int ret = cholmod_updown(UpOrDown, C_cs_perm, Base::m_cholmodFactor, &this->cholmod());
111  cholmod_free_sparse(&C_cs_perm, &this->cholmod());
112  assert(ret != 0);
113  return ret;
114  }
115 
117 
118 protected:
119  void init() {
120  m_cholmod.final_asis = 1;
121  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
122  // In CholmodBase::CholmodBase(), the following statement is missing in
123  // SuiteSparse 3.2.0.8. Fixed in 3.2.7
124  Base::m_shiftOffset[0] = Base::m_shiftOffset[1] = RealScalar(0.0);
125  }
126 };
127 
128 using namespace std;
129 
130 static double sqr(double x) { return x * x; }
131 
132 // const double posErrorIncrement=0.02;
133 
134 namespace {
135 LOG_LOGGER _log = LOG_GET("jointcal.AstrometryFit");
136 }
137 
138 namespace lsst {
139 namespace jointcal {
140 
141 AstrometryFit::AstrometryFit(Associations &associations, AstrometryModel *astrometryModel, double posError)
142  : _assoc(associations), _astrometryModel(astrometryModel), _posError(posError) {
143  _LastNTrip = 0;
144  _JDRef = 0;
145 
146  _posError = posError;
147 
148  _referenceColor = 0;
149  _sigCol = 0;
150  unsigned count = 0;
151  for (auto const &i : _assoc.fittedStarList) {
152  _referenceColor += i->color;
153  _sigCol += sqr(i->color);
154  count++;
155  }
156  if (count) {
157  _referenceColor /= double(count);
158  if (_sigCol > 0) _sigCol = sqrt(_sigCol / count - sqr(_referenceColor));
159  }
160  LOGLS_INFO(_log, "Reference Color: " << _referenceColor << " sig " << _sigCol);
161 
162  _nRefrac = _assoc.NBands();
163  _refractionCoefficient = 0;
164 
165  _nMeasuredStars = 0;
166  // The various _npar... are initialized in assignIndices.
167  // Although there is no reason to adress them before one might be tempted by
168  // evaluating a Chi2 rightaway, .. which uses these counts, so:
169  assignIndices("");
170 }
171 
172 #define NPAR_PM 2
173 
174 /* ! this routine is used in 3 instances: when computing
175 the derivatives, when computing the Chi2, when filling a tuple.
176 */
177 Point AstrometryFit::transformFittedStar(const FittedStar &fittedStar, const Gtransfo *sky2TP,
178  const Point &refractionVector, double refractionCoeff,
179  double mjd) const {
180  Point fittedStarInTP = sky2TP->apply(fittedStar);
181  if (fittedStar.mightMove) {
182  fittedStarInTP.x += fittedStar.pmx * mjd;
183  fittedStarInTP.y += fittedStar.pmy * mjd;
184  }
185  // account for atmospheric refraction: does nothing if color
186  // have not been assigned
187  // the color definition shouldbe the same when computing derivatives
188  double color = fittedStar.color - _referenceColor;
189  fittedStarInTP.x += refractionVector.x * color * refractionCoeff;
190  fittedStarInTP.y += refractionVector.y * color * refractionCoeff;
191  return fittedStarInTP;
192 }
193 
197 static void tweakAstromMeasurementErrors(FatPoint &P, const MeasuredStar &Ms, double error) {
198  static bool called = false;
199  static double increment = 0;
200  if (!called) {
201  increment = sqr(error); // was in Preferences
202  called = true;
203  }
204  P.vx += increment;
205  P.vy += increment;
206 }
207 
208 // we could consider computing the chi2 here.
209 // (although it is not extremely useful)
210 void AstrometryFit::LSDerivatives1(const CcdImage &ccdImage, TripletList &tList, Eigen::VectorXd &rhs,
211  const MeasuredStarList *msList) const {
212  /***************************************************************************/
214  /***************************************************************************/
215  /* Setup */
216  /* this routine works in two different ways: either providing the
217  ccdImage, of providing the MeasuredStarList. In the latter case, the
218  ccdImage should match the one(s) in the list. */
219  if (msList) assert(&(msList->front()->getCcdImage()) == &ccdImage);
220 
221  // get the Mapping
222  const Mapping *mapping = _astrometryModel->getMapping(ccdImage);
223  // count parameters
224  unsigned npar_mapping = (_fittingDistortions) ? mapping->getNpar() : 0;
225  unsigned npar_pos = (_fittingPos) ? 2 : 0;
226  unsigned npar_refrac = (_fittingRefrac) ? 1 : 0;
227  unsigned npar_pm = (_fittingPM) ? NPAR_PM : 0;
228  unsigned npar_tot = npar_mapping + npar_pos + npar_refrac + npar_pm;
229  // if (npar_tot == 0) this CcdImage does not contribute
230  // any constraint to the fit, so :
231  if (npar_tot == 0) return;
232  vector<unsigned> indices(npar_tot, -1);
233  if (_fittingDistortions) mapping->setMappingIndices(indices);
234 
235  // proper motion stuff
236  double mjd = ccdImage.getMjd() - _JDRef;
237  // refraction stuff
238  Point refractionVector = ccdImage.getRefractionVector();
239  // transformation from sky to TP
240  const Gtransfo *sky2TP = _astrometryModel->getSky2TP(ccdImage);
241  // reserve matrices once for all measurements
242  GtransfoLin dypdy;
243  // the shape of h (et al) is required this way in order to be able to
244  // separate derivatives along x and y as vectors.
245  Eigen::MatrixX2d h(npar_tot, 2), halpha(npar_tot, 2), hw(npar_tot, 2);
246  Eigen::Matrix2d transW(2, 2);
247  Eigen::Matrix2d alpha(2, 2);
248  Eigen::VectorXd grad(npar_tot);
249  // current position in the Jacobian
250  unsigned kTriplets = tList.getNextFreeIndex();
251  const MeasuredStarList &catalog = (msList) ? *msList : ccdImage.getCatalogForFit();
252 
253  for (auto &i : catalog) {
254  const MeasuredStar &ms = *i;
255  if (!ms.isValid()) continue;
256  // tweak the measurement errors
257  FatPoint inPos = ms;
258  tweakAstromMeasurementErrors(inPos, ms, _posError);
259  h.setZero(); // we cannot be sure that all entries will be overwritten.
260  FatPoint outPos;
261  // should *not* fill h if whatToFit excludes mapping parameters.
262  if (_fittingDistortions)
263  mapping->computeTransformAndDerivatives(inPos, outPos, h);
264  else
265  mapping->transformPosAndErrors(inPos, outPos);
266 
267  unsigned ipar = npar_mapping;
268  double det = outPos.vx * outPos.vy - sqr(outPos.vxy);
269  if (det <= 0 || outPos.vx <= 0 || outPos.vy <= 0) {
270  LOGLS_WARN(_log, "Inconsistent measurement errors: drop measurement at "
271  << Point(ms) << " in image " << ccdImage.getName());
272  continue;
273  }
274  transW(0, 0) = outPos.vy / det;
275  transW(1, 1) = outPos.vx / det;
276  transW(0, 1) = transW(1, 0) = -outPos.vxy / det;
277  // compute alpha, a triangular square root
278  // of transW (i.e. a Cholesky factor)
279  alpha(0, 0) = sqrt(transW(0, 0));
280  // checked that alpha*alphaT = transW
281  alpha(1, 0) = transW(0, 1) / alpha(0, 0);
282  // DB - I think that the next line is equivalent to : alpha(1,1) = 1./sqrt(outPos.vy)
283  // PA - seems correct !
284  alpha(1, 1) = 1. / sqrt(det * transW(0, 0));
285  alpha(0, 1) = 0;
286 
287  auto fs = ms.getFittedStar();
288 
289  Point fittedStarInTP =
290  transformFittedStar(*fs, sky2TP, refractionVector, _refractionCoefficient, mjd);
291 
292  // compute derivative of TP position w.r.t sky position ....
293  if (npar_pos > 0) // ... if actually fitting FittedStar position
294  {
295  sky2TP->computeDerivative(*fs, dypdy, 1e-3);
296  // sign checked
297  // TODO Still have to check with non trivial non-diagonal terms
298  h(npar_mapping, 0) = -dypdy.A11();
299  h(npar_mapping + 1, 0) = -dypdy.A12();
300  h(npar_mapping, 1) = -dypdy.A21();
301  h(npar_mapping + 1, 1) = -dypdy.A22();
302  indices[npar_mapping] = fs->getIndexInMatrix();
303  indices.at(npar_mapping + 1) = fs->getIndexInMatrix() + 1;
304  ipar += npar_pos;
305  }
306  /* only consider proper motions of objects allowed to move,
307  unless the fit is going to be degenerate */
308  if (_fittingPM && fs->mightMove) {
309  h(ipar, 0) = -mjd; // Sign unchecked but consistent with above
310  h(ipar + 1, 1) = -mjd;
311  indices[ipar] = fs->getIndexInMatrix() + 2;
312  indices[ipar + 1] = fs->getIndexInMatrix() + 3;
313  ipar += npar_pm;
314  }
315  if (_fittingRefrac) {
316  /* if the definition of color changes, it has to remain
317  consistent with transformFittedStar */
318  double color = fs->color - _referenceColor;
319  // sign checked
320  h(ipar, 0) = -refractionVector.x * color;
321  h(ipar, 1) = -refractionVector.y * color;
322  indices[ipar] = _refracPosInMatrix;
323  ipar += 1;
324  }
325 
326  // We can now compute the residual
327  Eigen::Vector2d res(fittedStarInTP.x - outPos.x, fittedStarInTP.y - outPos.y);
328 
329  // do not write grad = h*transW*res to avoid
330  // dynamic allocation of a temporary
331  halpha = h * alpha;
332  hw = h * transW;
333  grad = hw * res;
334  // now feed in triplets and rhs
335  for (unsigned ipar = 0; ipar < npar_tot; ++ipar) {
336  for (unsigned ic = 0; ic < 2; ++ic) {
337  double val = halpha(ipar, ic);
338  if (val == 0) continue;
339 #if (TRIPLET_INTERNAL_COORD == COL)
340  tList.addTriplet(indices[ipar], kTriplets + ic, val);
341 #else
342  tList.addTriplet(kTriplets + ic, indices[ipar], val);
343 #endif
344  }
345  rhs(indices[ipar]) += grad(ipar);
346  }
347  kTriplets += 2; // each measurement contributes 2 columns in the Jacobian
348  } // end loop on measurements
349  tList.setNextFreeIndex(kTriplets);
350 }
351 
352 // we could consider computing the chi2 here.
353 // (although it is not extremely useful)
354 
355 #define HACK_REF_ERRORS 1. // used to isolate the measurement or ref terms
356 
358  Eigen::VectorXd &rhs) const {
359  /* We compute here the derivatives of the terms involving fitted
360  stars and reference stars. They only provide contributions if we
361  are fitting positions: */
362  if (!_fittingPos) return;
363  /* the other case where the accumulation of derivatives stops
364  here is when there are no RefStars */
365  if (_assoc.refStarList.size() == 0) return;
366  Eigen::Matrix2d w(2, 2);
367  Eigen::Matrix2d alpha(2, 2);
368  Eigen::Matrix2d h(2, 2), halpha(2, 2), hw(2, 2);
369  GtransfoLin der;
370  Eigen::Vector2d res, grad;
371  unsigned indices[2 + NPAR_PM];
372  unsigned kTriplets = tList.getNextFreeIndex();
373  /* We cannot use the spherical coordinates directly to evaluate
374  Euclidean distances, we have to use a projector on some plane in
375  order to express least squares. Not projecting could lead to a
376  disaster around the poles or across alpha=0. So we need a
377  projector. We construct a projector and will change its
378  projection point at every object */
379  TanRaDec2Pix proj(GtransfoLin(), Point(0., 0.));
380  for (auto const &i : fsl) {
381  const FittedStar &fs = *i;
382  const RefStar *rs = fs.getRefStar();
383  if (rs == nullptr) continue;
384  proj.setTangentPoint(fs);
385  // fs projects to (0,0), no need to compute its transform.
386  FatPoint rsProj;
387  proj.transformPosAndErrors(*rs, rsProj);
388  proj.computeDerivative(fs, der, 1e-4);
389  // sign checked. TODO check that the off-diagonal terms are OK.
390  h(0, 0) = -der.A11();
391  h(1, 0) = -der.A12();
392  h(0, 1) = -der.A21();
393  h(1, 1) = -der.A22();
394  // TO DO : account for proper motions.
395  double det = rsProj.vx * rsProj.vy - sqr(rsProj.vxy);
396  if (rsProj.vx <= 0 || rsProj.vy <= 0 || det <= 0) {
397  LOGLS_WARN(_log, "RefStar error matrix not positive definite for: " << *rs);
398  continue;
399  }
400  w(0, 0) = rsProj.vy / det;
401  w(0, 1) = w(1, 0) = -rsProj.vxy / det;
402  w(1, 1) = rsProj.vx / det;
403  w *= HACK_REF_ERRORS;
404  det /= sqr(HACK_REF_ERRORS);
405  // compute alpha, a triangular square root
406  // of w (i.e. a Cholesky factor)
407  alpha(0, 0) = sqrt(w(0, 0));
408  // checked that alpha*alphaT = transW
409  alpha(1, 0) = w(0, 1) / alpha(0, 0);
410  alpha(1, 1) = 1. / sqrt(det * w(0, 0));
411  alpha(0, 1) = 0;
412  indices[0] = fs.getIndexInMatrix();
413  indices[1] = fs.getIndexInMatrix() + 1;
414  unsigned npar_tot = 2;
415  /* TODO: account here for proper motions in the reference
416  catalog. We can code the effect and set the value to 0. Most
417  (all?) catalogs do not even come with a reference epoch. Gaia
418  will change that. When refraction enters into the game, one should
419  pay attention to the orientation of the frame */
420 
421  /* The residual should be Proj(fs)-Proj(*rs) in order to be consistent
422  with the measurement terms. Since P(fs) = 0, we have: */
423  res[0] = -rsProj.x;
424  res[1] = -rsProj.y;
425  halpha = h * alpha;
426  // grad = h*w*res
427  hw = h * w;
428  grad = hw * res;
429  // now feed in triplets and rhs
430  for (unsigned ipar = 0; ipar < npar_tot; ++ipar) {
431  for (unsigned ic = 0; ic < 2; ++ic) {
432  double val = halpha(ipar, ic);
433  if (val == 0) continue;
434 #if (TRIPLET_INTERNAL_COORD == COL)
435  tList.addTriplet(indices[ipar], kTriplets + ic, val);
436 #else
437  tList.addTriplet(kTriplets + ic, indices[ipar], val);
438 #endif
439  }
440  rhs(indices[ipar]) += grad(ipar);
441  }
442  kTriplets += 2; // each measurement contributes 2 columns in the Jacobian
443  }
444  tList.setNextFreeIndex(kTriplets);
445 }
446 
449 void AstrometryFit::LSDerivatives(TripletList &tList, Eigen::VectorXd &rhs) const {
450  auto L = _assoc.getCcdImageList();
451  for (auto const &im : L) {
452  LSDerivatives1(*im, tList, rhs);
453  }
454  LSDerivatives2(_assoc.fittedStarList, tList, rhs);
455 }
456 
457 // This is almost a selection of lines of LSDerivatives1(CcdImage ...)
458 /* This routine (and the following one) is template because it is used
459 both with its first argument as "const CCdImage &" and "CcdImage &",
460 and I did not want to replicate it. The constness of the iterators is
461 automagically set by declaring them as "auto" */
462 
463 template <class ImType, class Accum>
464 void AstrometryFit::accumulateStatImage(ImType &image, Accum &accu) const {
465  /**********************************************************************/
467  /**********************************************************************/
468  /* Setup */
469  // 1 : get the Mapping's
470  const Mapping *mapping = _astrometryModel->getMapping(image);
471  // proper motion stuff
472  double mjd = image.getMjd() - _JDRef;
473  // refraction stuff
474  Point refractionVector = image.getRefractionVector();
475  // transformation from sky to TP
476  const Gtransfo *sky2TP = _astrometryModel->getSky2TP(image);
477  // reserve matrix once for all measurements
478  Eigen::Matrix2Xd transW(2, 2);
479 
480  auto &catalog = image.getCatalogForFit();
481  for (auto const &ms : catalog) {
482  if (!ms->isValid()) continue;
483  // tweak the measurement errors
484  FatPoint inPos = *ms;
485  tweakAstromMeasurementErrors(inPos, *ms, _posError);
486 
487  FatPoint outPos;
488  // should *not* fill h if whatToFit excludes mapping parameters.
489  mapping->transformPosAndErrors(inPos, outPos);
490  double det = outPos.vx * outPos.vy - sqr(outPos.vxy);
491  if (det <= 0 || outPos.vx <= 0 || outPos.vy <= 0) {
492  LOGLS_WARN(_log, " Inconsistent measurement errors :drop measurement at "
493  << Point(*ms) << " in image " << image.getName());
494  continue;
495  }
496  transW(0, 0) = outPos.vy / det;
497  transW(1, 1) = outPos.vx / det;
498  transW(0, 1) = transW(1, 0) = -outPos.vxy / det;
499 
500  auto fs = ms->getFittedStar();
501  Point fittedStarInTP =
502  transformFittedStar(*fs, sky2TP, refractionVector, _refractionCoefficient, mjd);
503 
504  Eigen::Vector2d res(fittedStarInTP.x - outPos.x, fittedStarInTP.y - outPos.y);
505  double chi2Val = res.transpose() * transW * res;
506 
507  accu.addEntry(chi2Val, 2, ms);
508  } // end of loop on measurements
509 }
510 
512 template <class ListType, class Accum>
513 void AstrometryFit::accumulateStatImageList(ListType &list, Accum &accum) const {
514  for (auto &im : list) {
515  accumulateStatImage(*im, accum);
516  }
517 }
518 
519 template <class Accum>
520 void AstrometryFit::accumulateStatRefStars(Accum &accum) const {
521  /* If you wonder why we project here, read comments in
522  AstrometryFit::LSDerivatives2(TripletList &TList, Eigen::VectorXd &Rhs) */
523  FittedStarList &fsl = _assoc.fittedStarList;
524  TanRaDec2Pix proj(GtransfoLin(), Point(0., 0.));
525  for (auto const &fs : fsl) {
526  const RefStar *rs = fs->getRefStar();
527  if (rs == nullptr) continue;
528  proj.setTangentPoint(*fs);
529  // fs projects to (0,0), no need to compute its transform.
530  FatPoint rsProj;
531  proj.transformPosAndErrors(*rs, rsProj);
532  // TO DO : account for proper motions.
533  double rx = rsProj.x; // -fsProj.x (which is 0)
534  double ry = rsProj.y;
535  double det = rsProj.vx * rsProj.vy - sqr(rsProj.vxy);
536  double wxx = rsProj.vy / det;
537  double wyy = rsProj.vx / det;
538  double wxy = -rsProj.vxy / det;
539  wxx *= HACK_REF_ERRORS;
540  wyy *= HACK_REF_ERRORS;
541  wxy *= HACK_REF_ERRORS;
542  accum.addEntry(wxx * sqr(rx) + 2 * wxy * rx * ry + wyy * sqr(ry), 2, fs);
543  }
544 }
545 
547  Chi2 chi2;
548  accumulateStatImageList(_assoc.getCcdImageList(), chi2);
549  // now ref stars:
550  accumulateStatRefStars(chi2);
551  // so far, ndof contains the number of squares.
552  // So, subtract here the number of parameters.
553  chi2.ndof -= _nParTot;
554  return chi2;
555 }
556 
558 
562 struct Chi2Entry {
563  double chi2;
564  std::shared_ptr<BaseStar> ps;
565 
566  Chi2Entry(double c, std::shared_ptr<BaseStar> s) : chi2(c), ps(std::move(s)) {}
567  // for sort
568  bool operator<(const Chi2Entry &R) const { return (chi2 < R.chi2); }
569 };
570 
571 struct Chi2Vect : public std::vector<Chi2Entry> {
572  void addEntry(double Chi2Val, unsigned ndof, std::shared_ptr<BaseStar> ps) {
573  this->push_back(Chi2Entry(Chi2Val, std::move(ps)));
574  }
575 };
576 
578 
580 void AstrometryFit::setMeasuredStarIndices(const MeasuredStar &ms, std::vector<unsigned> &indices) const {
581  if (_fittingDistortions) {
582  const Mapping *mapping = _astrometryModel->getMapping(ms.getCcdImage());
583  mapping->setMappingIndices(indices);
584  }
585  auto fs = ms.getFittedStar();
586  unsigned fsIndex = fs->getIndexInMatrix();
587  if (_fittingPos) {
588  indices.push_back(fsIndex);
589  indices.push_back(fsIndex + 1);
590  }
591  // For securing the outlier removal, the next block is just useless
592  if (_fittingPM) {
593  for (unsigned k = 0; k < NPAR_PM; ++k) indices.push_back(fsIndex + 2 + k);
594  }
595  /* Should not put the index of refaction stuff or we will not be
596  able to remove more than 1 star at a time. */
597 }
598 
601  TripletList &tList, Eigen::VectorXd &grad) {
602  // contributions from measurement terms:
603  for (auto const &i : msOutliers) {
604  MeasuredStarList tmp;
605  tmp.push_back(i);
606  const CcdImage &ccd = i->getCcdImage();
607  LSDerivatives1(ccd, tList, grad, &tmp);
608  }
609  LSDerivatives2(fOutliers, tList, grad);
610 }
611 
612 unsigned AstrometryFit::removeOutliers(double nSigmaCut, const std::string &measOrRef) {
613  MeasuredStarList msOutliers;
614  FittedStarList fsOutliers;
615  unsigned n = findOutliers(nSigmaCut, msOutliers, fsOutliers, measOrRef);
616  removeMeasOutliers(msOutliers);
617  removeRefOutliers(fsOutliers);
618  return n;
619 }
620 
621 unsigned AstrometryFit::findOutliers(double nSigmaCut, MeasuredStarList &msOutliers,
622  FittedStarList &fsOutliers, const std::string &measOrRef) const {
623  bool searchMeas = (measOrRef.find("Meas") != std::string::npos);
624  bool searchRef = (measOrRef.find("Ref") != std::string::npos);
625 
626  // collect chi2 contributions
627  Chi2Vect chi2s;
628  chi2s.reserve(_nMeasuredStars + _assoc.refStarList.size());
629  // contributions from measurement terms:
630  if (searchMeas) accumulateStatImageList(_assoc.ccdImageList, chi2s);
631  // and from reference terms
632  if (searchRef) accumulateStatRefStars(chi2s);
633 
634  // do some stat
635  unsigned nval = chi2s.size();
636  if (nval == 0) return 0;
637  sort(chi2s.begin(), chi2s.end());
638  double median =
639  (nval & 1) ? chi2s[nval / 2].chi2 : 0.5 * (chi2s[nval / 2 - 1].chi2 + chi2s[nval / 2].chi2);
640  // some more stats. should go into the class if recycled anywhere else
641  double sum = 0;
642  double sum2 = 0;
643  for (auto i = chi2s.begin(); i != chi2s.end(); ++i) {
644  sum += i->chi2;
645  sum2 += sqr(i->chi2);
646  }
647  double average = sum / nval;
648  double sigma = sqrt(sum2 / nval - sqr(average));
649  LOGLS_DEBUG(_log,
650  "RemoveOutliers chi2 stat: mean/median/sigma " << average << '/' << median << '/' << sigma);
651  double cut = average + nSigmaCut * sigma;
652  /* For each of the parameters, we will not remove more than 1
653  measurement that contributes to constraining it. Keep track using
654  of what we are touching using an integer vector. This is the
655  trick that Marc Betoule came up to for outlier removals in "star
656  flats" fits. */
657  Eigen::VectorXi affectedParams(_nParTot);
658  affectedParams.setZero();
659 
660  unsigned nOutliers = 0; // returned to the caller
661  // start from the strongest outliers.
662  for (auto i = chi2s.rbegin(); i != chi2s.rend(); ++i) {
663  if (i->chi2 < cut) break; // because the array is sorted.
664  vector<unsigned> indices;
665  indices.reserve(100); // just there to limit reallocations.
666  /* now, we want to get the indices of the parameters this chi2
667  term depends on. We have to figure out which kind of term it
668  is; we use for that the type of the star attached to
669  the Chi2Entry. */
670  auto ms = std::dynamic_pointer_cast<MeasuredStar>(i->ps);
671  std::shared_ptr<FittedStar> fs;
672  if (!ms) // it is reference term.
673  {
674  fs = std::dynamic_pointer_cast<FittedStar>(i->ps);
675  indices.push_back(fs->getIndexInMatrix());
676  indices.push_back(fs->getIndexInMatrix() + 1); // probably useless
677  /* One might think it would be useful to account for PM
678  parameters here, but it is just useless */
679  } else // it is a measurement term.
680  {
681  setMeasuredStarIndices(*ms, indices);
682  }
683 
684  /* Find out if we already discarded a stronger outlier
685  constraining some parameter this one constrains as well. If
686  yes, we keep this one, because this stronger outlier could be
687  causing the large chi2 we have in hand. */
688  bool drop_it = true;
689  for (auto const &i : indices)
690  if (affectedParams(i) != 0) drop_it = false;
691 
692  if (drop_it) // store the outlier in one of the lists:
693  {
694  if (ms) // measurement term
695  msOutliers.push_back(ms);
696  else // ref term
697  fsOutliers.push_back(fs);
698  // mark the parameters as directly changed when we discard this chi2 term.
699  for (auto const &i : indices) affectedParams(i)++;
700  nOutliers++;
701  }
702  } // end loop on measurements/references
703  LOGLS_INFO(_log, "findOutliers: found " << msOutliers.size() << " meas outliers and " << fsOutliers.size()
704  << " ref outliers ");
705 
706  return nOutliers;
707 }
708 
710  for (auto &i : outliers) {
711  MeasuredStar &ms = *i;
712  auto fs = std::const_pointer_cast<FittedStar>(ms.getFittedStar());
713  ms.setValid(false);
714  fs->getMeasurementCount()--; // could be put in setValid
715  }
716 }
717 
719  for (auto &i : outliers) {
720  FittedStar &fs = *i;
721  fs.setRefStar(nullptr);
722  }
723 }
724 
725 void AstrometryFit::assignIndices(const std::string &whatToFit) {
726  _whatToFit = whatToFit;
727  LOGLS_INFO(_log, "assignIndices: Now fitting " << whatToFit);
728  _fittingDistortions = (_whatToFit.find("Distortions") != string::npos);
729  _fittingPos = (_whatToFit.find("Positions") != string::npos);
730  _fittingRefrac = (_whatToFit.find("Refrac") != string::npos);
731  if (_sigCol == 0 && _fittingRefrac) {
732  LOGLS_WARN(_log,
733  "Cannot fit refraction coefficients without a color lever arm. Ignoring refraction.");
734  _fittingRefrac = false;
735  }
736  _fittingPM = (_whatToFit.find("PM") != string::npos);
737  // When entering here, we assume that whatToFit has already been interpreted.
738 
739  _nParDistortions = 0;
740  if (_fittingDistortions) _nParDistortions = _astrometryModel->assignIndices(0, _whatToFit);
741  unsigned ipar = _nParDistortions;
742 
743  if (_fittingPos) {
744  FittedStarList &fsl = _assoc.fittedStarList;
745  for (auto const &i : fsl) {
746  FittedStar &fs = *i;
747  // the parameter layout here is used also
748  // - when filling the derivatives
749  // - when updating (offsetParams())
750  // - in GetMeasuredStarIndices
751  fs.setIndexInMatrix(ipar);
752  ipar += 2;
753  if ((_fittingPM)&fs.mightMove) ipar += NPAR_PM;
754  }
755  }
756  _nParPositions = ipar - _nParDistortions;
757  if (_fittingRefrac) {
758  _refracPosInMatrix = ipar;
759  ipar += _nRefrac;
760  }
761  _nParTot = ipar;
762 }
763 
764 void AstrometryFit::offsetParams(const Eigen::VectorXd &delta) {
765  if (delta.size() != _nParTot)
766  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
767  "AstrometryFit::offsetParams : the provided vector length is not compatible with "
768  "the current whatToFit setting");
769  if (_fittingDistortions) _astrometryModel->offsetParams(delta);
770 
771  if (_fittingPos) {
772  FittedStarList &fsl = _assoc.fittedStarList;
773  for (auto const &i : fsl) {
774  FittedStar &fs = *i;
775  // the parameter layout here is used also
776  // - when filling the derivatives
777  // - when assigning indices (assignIndices())
778  unsigned index = fs.getIndexInMatrix();
779  fs.x += delta(index);
780  fs.y += delta(index + 1);
781  if ((_fittingPM)&fs.mightMove) {
782  fs.pmx += delta(index + 2);
783  fs.pmy += delta(index + 3);
784  }
785  }
786  }
787  if (_fittingRefrac) {
788  _refractionCoefficient += delta(_refracPosInMatrix);
789  }
790 }
791 
792 // should not be too large !
793 #ifdef STORAGE
794 static void write_sparse_matrix_in_fits(const SpMat &mat, const string &fitsName) {
795  if (mat.rows() * mat.cols() > 2e8) {
796  LOGLS_WARN(_log,
797  "write_sparse_matrix_in_fits: yout matrix is too large. " << fitsName << " not generated");
798  return;
799  }
800  Mat m(mat.rows(), mat.cols());
801  for (int k = 0; k < mat.outerSize(); ++k)
802  for (SpMat::InnerIterator it(mat, k); it; ++it) {
803  m(it.row(), it.col()) = it.value();
804  }
805  m.writeFits(fitsName);
806 }
807 
808 static void write_vect_in_fits(const Eigen::VectorXd &vectorXd, const string &fitsName) {
809  Vect v(vectorXd.size());
810  for (int k = 0; k < vectorXd.size(); ++k) v(k) = V(k);
811  Mat(v).writeFits(fitsName);
812 }
813 
814 #endif
815 
816 unsigned AstrometryFit::minimize(const std::string &whatToFit, const double nSigRejCut) {
817  assignIndices(whatToFit);
818 
819  // return code can take 3 values :
820  // 0 : fit has converged - no more outliers
821  // 1 : still some ouliers but chi2 increases
822  // 2 : factorization failed
823  unsigned returnCode = 0;
824 
825  // TODO : write a guesser for the number of triplets
826  unsigned nTrip = (_LastNTrip) ? _LastNTrip : 1e6;
827  TripletList tList(nTrip);
828  Eigen::VectorXd grad(_nParTot);
829  grad.setZero();
830 
831  // Fill the triplets
832  LSDerivatives(tList, grad);
833  _LastNTrip = tList.size();
834 
835  LOGLS_DEBUG(_log, "End of triplet filling, ntrip = " << tList.size());
836 
837  SpMat hessian;
838  {
839 #if (TRIPLET_INTERNAL_COORD == COL)
840  SpMat jacobian(_nParTot, tList.getNextFreeIndex());
841  jacobian.setFromTriplets(tList.begin(), tList.end());
842  // release memory shrink_to_fit is C++11
843  tList.clear(); // tList.shrink_to_fit();
844  hessian = jacobian * jacobian.transpose();
845 #else
846  SpMat jacobian(tList.NextRank(), _nParTot);
847  jacobian.setFromTriplets(tList.begin(), tList.end());
848  // release memory shrink_to_fit is C++11
849  tList.clear(); // tList.shrink_to_fit();
850  LOGLS_DEBUG(_log, " starting H=JtJ ");
851  hessian = jacobian.transpose() * jacobian;
852 #endif
853  } // release the Jacobian
854 
855  LOGLS_DEBUG(_log, "Starting factorization, hessian: dim="
856  << hessian.rows() << " nnz=" << hessian.nonZeros()
857  << " filling-frac = " << hessian.nonZeros() / sqr(hessian.rows()));
858 
859  CholmodSimplicialLDLT2<SpMat> chol(hessian);
860  if (chol.info() != Eigen::Success) {
861  LOGLS_ERROR(_log, "minimize: factorization failed ");
862  return 2;
863  }
864 
865  unsigned tot_outliers = 0;
866  double old_chi2 = computeChi2().chi2;
867 
868  while (true) {
869  Eigen::VectorXd delta = chol.solve(grad);
870  offsetParams(delta);
871  Chi2 current_chi2(computeChi2());
872  LOGLS_DEBUG(_log, current_chi2);
873  if (current_chi2.chi2 > old_chi2) {
874  LOGL_WARN(_log, "chi2 went up, exiting outlier rejection loop");
875  returnCode = 1;
876  break;
877  }
878  old_chi2 = current_chi2.chi2;
879 
880  if (nSigRejCut == 0) break;
881  MeasuredStarList moutliers;
882  FittedStarList foutliers;
883  int n_outliers = findOutliers(nSigRejCut, moutliers, foutliers);
884  tot_outliers += n_outliers;
885  if (n_outliers == 0) break;
886  TripletList tList(1000); // initial allocation size.
887  grad.setZero(); // recycle the gradient
888  // compute the contributions of outliers to derivatives
889  outliersContributions(moutliers, foutliers, tList, grad);
890  // actually discard them
891  removeMeasOutliers(moutliers);
892  removeRefOutliers(foutliers);
893  // convert triplet list to eigen internal format
894  SpMat h(_nParTot, tList.getNextFreeIndex());
895  h.setFromTriplets(tList.begin(), tList.end());
896  int update_status = chol.update(h, false /* means downdate */);
897  LOGLS_DEBUG(_log, "cholmod update_status " << update_status);
898  /* The contribution of outliers to the gradient is the opposite
899  of the contribution of all other terms, because they add up
900  to 0 */
901  grad *= -1;
902  }
903 
904  LOGLS_INFO(_log, "Total number of outliers " << tot_outliers);
905 
906  return returnCode;
907 }
908 
910 #if (0)
911  const char *what2fit[] = {"Positions",
912  "Distortions",
913  "Refrac",
914  "Positions Distortions",
915  "Positions Refrac",
916  "Distortions Refrac",
917  "Positions Distortions Refrac"};
918 #endif
919  const char *what2fit[] = {"Positions", "Distortions", "Positions Distortions"};
920  // DEBUG
921  for (unsigned k = 0; k < sizeof(what2fit) / sizeof(what2fit[0]); ++k) {
922  assignIndices(what2fit[k]);
923  TripletList tList(10000);
924  Eigen::VectorXd rhs(_nParTot);
925  rhs.setZero();
926  LSDerivatives(tList, rhs);
927  SpMat jacobian(_nParTot, tList.getNextFreeIndex());
928  jacobian.setFromTriplets(tList.begin(), tList.end());
929  SpMat hessian = jacobian * jacobian.transpose();
930 #ifdef STORAGE
931  char name[24];
932  sprintf(name, "h%d.fits", k);
933  write_sparse_matrix_in_fits(hessian, name);
934  sprintf(name, "g%d.fits", k);
935  write_vect_in_fits(rhs, name);
936 #endif
937  LOGLS_DEBUG(_log, "npar : " << _nParTot << ' ' << _nParDistortions);
938  }
939 }
940 
941 void AstrometryFit::makeResTuple(const std::string &tupleName) const {
942  /* cook-up 2 different file names by inserting something just before
943  the dot (if any), and within the actual file name. */
944  size_t dot = tupleName.rfind('.');
945  size_t slash = tupleName.rfind('/');
946  if (dot == string::npos || (slash != string::npos && dot < slash)) dot = tupleName.size();
947  std::string meas_tuple(tupleName);
948  meas_tuple.insert(dot, "-meas");
949  makeMeasResTuple(meas_tuple);
950  std::string ref_tuple(tupleName);
951  ref_tuple.insert(dot, "-ref");
952  makeRefResTuple(ref_tuple);
953 }
954 
955 void AstrometryFit::makeMeasResTuple(const std::string &tupleName) const {
956  std::ofstream tuple(tupleName.c_str());
957  tuple << "#xccd: coordinate in CCD" << endl
958  << "#yccd: " << endl
959  << "#rx: residual in degrees in TP" << endl
960  << "#ry:" << endl
961  << "#xtp: transformed coordinate in TP " << endl
962  << "#ytp:" << endl
963  << "#mag: rough mag" << endl
964  << "#jd: Julian date of the measurement" << endl
965  << "#rvx: transformed measurement uncertainty " << endl
966  << "#rvy:" << endl
967  << "#rvxy:" << endl
968  << "#color : " << endl
969  << "#fsindex: some unique index of the object" << endl
970  << "#ra: pos of fitted star" << endl
971  << "#dec: pos of fitted star" << endl
972  << "#chi2: contribution to Chi2 (2D dofs)" << endl
973  << "#nm: number of measurements of this FittedStar" << endl
974  << "#chip: chip number" << endl
975  << "#visit: visit id" << endl
976  << "#end" << endl;
977  const CcdImageList &L = _assoc.getCcdImageList();
978  for (auto const &i : L) {
979  const CcdImage &im = *i;
980  const MeasuredStarList &cat = im.getCatalogForFit();
981  const Mapping *mapping = _astrometryModel->getMapping(im);
982  const Point &refractionVector = im.getRefractionVector();
983  double mjd = im.getMjd() - _JDRef;
984  for (auto const &is : cat) {
985  const MeasuredStar &ms = *is;
986  if (!ms.isValid()) continue;
987  FatPoint tpPos;
988  FatPoint inPos = ms;
989  tweakAstromMeasurementErrors(inPos, ms, _posError);
990  mapping->transformPosAndErrors(inPos, tpPos);
991  const Gtransfo *sky2TP = _astrometryModel->getSky2TP(im);
992  auto fs = ms.getFittedStar();
993 
994  Point fittedStarInTP =
995  transformFittedStar(*fs, sky2TP, refractionVector, _refractionCoefficient, mjd);
996  Point res = tpPos - fittedStarInTP;
997  double det = tpPos.vx * tpPos.vy - sqr(tpPos.vxy);
998  double wxx = tpPos.vy / det;
999  double wyy = tpPos.vx / det;
1000  double wxy = -tpPos.vxy / det;
1001  // double chi2 = rx*(wxx*rx+wxy*ry)+ry*(wxy*rx+wyy*ry);
1002  double chi2 = wxx * res.x * res.x + wyy * res.y * res.y + 2 * wxy * res.x * res.y;
1003  tuple << std::setprecision(9);
1004  tuple << ms.x << ' ' << ms.y << ' ' << res.x << ' ' << res.y << ' ' << tpPos.x << ' ' << tpPos.y
1005  << ' ' << fs->getMag() << ' ' << mjd << ' ' << tpPos.vx << ' ' << tpPos.vy << ' '
1006  << tpPos.vxy << ' ' << fs->color << ' ' << fs->getIndexInMatrix() << ' ' << fs->x << ' '
1007  << fs->y << ' ' << chi2 << ' ' << fs->getMeasurementCount() << ' ' << im.getCcdId() << ' '
1008  << im.getVisit() << endl;
1009  } // loop on measurements in image
1010  } // loop on images
1011 }
1012 
1013 void AstrometryFit::makeRefResTuple(const std::string &tupleName) const {
1014  std::ofstream tuple(tupleName.c_str());
1015  tuple << "#ra: coordinates of FittedStar" << endl
1016  << "#dec: " << endl
1017  << "#rx: residual in degrees in TP" << endl
1018  << "#ry:" << endl
1019  << "#mag: mag" << endl
1020  << "#rvx: transformed measurement uncertainty " << endl
1021  << "#rvy:" << endl
1022  << "#rvxy:" << endl
1023  << "#color : " << endl
1024  << "#fsindex: some unique index of the object" << endl
1025  << "#chi2: contribution to Chi2 (2D dofs)" << endl
1026  << "#nm: number of measurements of this FittedStar" << endl
1027  << "#end" << endl;
1028  // The following loop is heavily inspired from AstrometryFit::computeChi2()
1029  const FittedStarList &fsl = _assoc.fittedStarList;
1030  TanRaDec2Pix proj(GtransfoLin(), Point(0., 0.));
1031  for (auto const &i : fsl) {
1032  const FittedStar &fs = *i;
1033  const RefStar *rs = fs.getRefStar();
1034  if (rs == nullptr) continue;
1035  proj.setTangentPoint(fs);
1036  // fs projects to (0,0), no need to compute its transform.
1037  FatPoint rsProj;
1038  proj.transformPosAndErrors(*rs, rsProj);
1039  double rx = rsProj.x; // -fsProj.x (which is 0)
1040  double ry = rsProj.y;
1041  double det = rsProj.vx * rsProj.vy - sqr(rsProj.vxy);
1042  double wxx = rsProj.vy / det;
1043  double wyy = rsProj.vx / det;
1044  double wxy = -rsProj.vxy / det;
1045  double chi2 = wxx * sqr(rx) + 2 * wxy * rx * ry + wyy * sqr(ry);
1046  tuple << std::setprecision(9);
1047  tuple << fs.x << ' ' << fs.y << ' ' << rx << ' ' << ry << ' ' << fs.getMag() << ' ' << rsProj.vx
1048  << ' ' << rsProj.vy << ' ' << rsProj.vxy << ' ' << fs.color << ' ' << fs.getIndexInMatrix()
1049  << ' ' << chi2 << ' ' << fs.getMeasurementCount() << endl;
1050  } // loop on FittedStars
1051 }
1052 } // namespace jointcal
1053 } // namespace lsst
Point getRefractionVector() const
Definition: CcdImage.h:159
void addEntry(double Chi2Val, unsigned ndof, std::shared_ptr< BaseStar > ps)
Objects used as position anchors, typically USNO stars.
Definition: RefStar.h:16
void outliersContributions(MeasuredStarList &mOutliers, FittedStarList &fOutLiers, TripletList &tList, Eigen::VectorXd &grad)
Contributions to derivatives from (presumably) outlier terms.
implements the linear transformations (6 real coefficients).
Definition: Gtransfo.h:288
void makeMeasResTuple(const std::string &tupleName) const
Produces a tuple containing residuals of measurement terms.
double getMjd() const
Julian Date.
Definition: CcdImage.h:136
virtual void transformPosAndErrors(const FatPoint &where, FatPoint &outPoint) const =0
The same as above but without the parameter derivatives (used to evaluate chi^2)
a class to accumulate chi2 contributions together with pointers to the contributors.
unsigned findOutliers(double nSigmaCut, MeasuredStarList &mSOutliers, FittedStarList &fSOutliers, const std::string &measOrRef="Meas Ref") const
Find Measurements and references contributing more than a cut, computed as The outliers are NOT remo...
virtual void setMappingIndices(std::vector< unsigned > &indices) const =0
Sets how this set of parameters (of length Npar()) map into the &quot;grand&quot; fit Expects that indices has ...
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:123
void makeResTuple(const std::string &tupleName) const
Produces both ntuples (cook up names from the provided string)
A point in a plane.
Definition: Point.h:13
void LSDerivatives(TripletList &tList, Eigen::VectorXd &rhs) const
Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting...
double A12() const
Definition: Gtransfo.h:327
std::shared_ptr< const FittedStar > getFittedStar() const
Definition: MeasuredStar.h:56
std::shared_ptr< BaseStar > ps
Interface class between AstrometryFit and an actual model for the Mapping (s) from pixels to some tan...
double chi2
unsigned removeOutliers(double nSigmaCut, const std::string &measOrRef="Meas Ref")
Discards measurements and reference contributions contributing to the chi2 more than nSigmaCut...
void removeMeasOutliers(MeasuredStarList &outliers)
Just removes outliers from the fit. No Refit done.
double getMag() const
derived using available zero points in input images. In the absence ofZP, ZP= 0.
Definition: FittedStar.h:102
The class that implements the relations between MeasuredStar and FittedStar.
Definition: Associations.h:28
A Point with uncertainties.
Definition: FatPoint.h:11
void LSDerivatives1(const CcdImage &ccdImage, TripletList &tList, Eigen::VectorXd &rhs, const MeasuredStarList *msList=nullptr) const
Compute derivatives of measurement terms for this CcdImage.
unsigned ndof
Definition: Chi2.h:15
unsigned minimize(const std::string &whatToFit, const double nSigRejCut=0)
Does a 1 step minimization, assuming a linear model.
A list of MeasuredStar. They are usually filled in Associations::AddImage.
Definition: MeasuredStar.h:71
double x
coordinate
Definition: Point.h:18
double A22() const
Definition: Gtransfo.h:329
Simple structure to accumulate Chi2 and Ndof.
Definition: Chi2.h:13
bool operator<(const Chi2Entry &R) const
virtual const Gtransfo * getSky2TP(const CcdImage &ccdImage) const =0
The transformation used to project the positions of FittedStars.
virtual const Mapping * getMapping(const CcdImage &) const =0
Mapping associated to a given CcdImage.
const CcdImageList & getCcdImageList() const
Definition: Associations.h:109
void setIndexInMatrix(const unsigned &index)
index is a value that a fit can set and reread....
Definition: FittedStar.h:114
void setRefStar(const RefStar *_refStar)
Set the astrometric reference star associated with this star.
Definition: FittedStar.cc:34
int getIndexInMatrix() const
Definition: FittedStar.h:117
void addTriplet(const unsigned i, const unsigned j, double val)
Definition: Tripletlist.h:23
void makeRefResTuple(const std::string &tupleName) const
Produces a tuple containing residuals of reference terms.
MatrixType::RealScalar RealScalar
CholmodSimplicialLDLT2(const MatrixType &matrix)
#define HACK_REF_ERRORS
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition: FittedStar.h:134
MatrixType::Index Index
Chi2Entry(double c, std::shared_ptr< BaseStar > s)
Chi2 computeChi2() const
Returns a chi2 for the current state.
double A21() const
Definition: Gtransfo.h:328
objects measured on actual images.
Definition: MeasuredStar.h:18
const MeasuredStarList & getCatalogForFit() const
Gets the catalog to be used for fitting, which may have been cleaned-up.
Definition: CcdImage.h:93
Eigen::Matrix< double, Eigen::Dynamic, 2 > MatrixX2d
Definition: Eigenstuff.h:6
unsigned NBands() const
Number of different bands in the input image list. Not implemented so far.
Definition: Associations.h:112
void transformPosAndErrors(const FatPoint &in, FatPoint &out) const
transform with analytical derivatives
Definition: Gtransfo.cc:1465
int update(const SpMat &H, const bool UpOrDown)
void assignIndices(const std::string &whatToFit)
Set parameter groups fixed or variable and assign indices to each parameter in the big matrix (which ...
virtual unsigned getNpar() const =0
Number of parameters in total.
#define NPAR_PM
void removeRefOutliers(FittedStarList &outliers)
Just removes outliers from the fit. No Refit done.
void setTangentPoint(const Point &tangentPoint)
Resets the projection (or tangent) point.
Definition: Gtransfo.cc:1445
unsigned getNextFreeIndex() const
Definition: Tripletlist.h:25
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane) ...
Definition: Gtransfo.h:513
Cholesky factorization class using cholmod, with the small-rank update capability.
void setValid(bool v)
Fits may use that to discard outliers.
Definition: MeasuredStar.h:65
a virtual (interface) class for geometric transformations.
Definition: Gtransfo.h:37
int getMeasurementCount() const
Definition: FittedStar.h:98
virtual void computeDerivative(const Point &where, GtransfoLin &derivative, const double step=0.01) const
Computes the local Derivative of a transfo. Step is used for numerical derivation.
Definition: Gtransfo.cc:75
std::list< std::shared_ptr< CcdImage > > CcdImageList
Definition: CcdImage.h:22
double A11() const
Definition: Gtransfo.h:326
Handler of an actual image from a single CCD.
Definition: CcdImage.h:31
Eigen::SparseMatrix< double > SpMat
virtual void computeTransformAndDerivatives(const FatPoint &where, FatPoint &outPoint, Eigen::MatrixX2d &H) const =0
Actually applies the mapping and evaluates the derivatives w.r.t the fitted parameters.
void checkStuff()
DEBUGGING routine.
bool isValid() const
Fits may use that to discard outliers.
Definition: MeasuredStar.h:63
void offsetParams(const Eigen::VectorXd &delta)
Offset the parameters by the requested quantities.
virtual unsigned assignIndices(unsigned firstIndex, const std::string &whatToFit)=0
Assign indices to parameters involved in mappings, starting at firstIndex.
std::string getName() const
Return the _name that identifies this ccdImage.
Definition: CcdImage.h:78
The objects which have been measured several times.
Definition: FittedStar.h:34
void setNextFreeIndex(unsigned index)
Definition: Tripletlist.h:27
VisitIdType getVisit() const
returns visit ID
Definition: CcdImage.h:130
virtual void offsetParams(const Eigen::VectorXd &delta)=0
Offset the parameters by the provided amounts.
int getCcdId() const
returns ccd ID
Definition: CcdImage.h:127
FittedStarList fittedStarList
Definition: Associations.h:32
void LSDerivatives2(const FittedStarList &fsl, TripletList &tList, Eigen::VectorXd &rhs) const
Compute derivatives of reference terms (if any), associated to the FittedStarList.