lsst.jointcal  master-gc935ebf72c
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
PhotometryFit.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 #include "Eigen/Sparse"
10 //#include "Eigen/CholmodSupport" // to switch to cholmod
11 #include "lsst/pex/exceptions.h"
12 #include <fstream>
14 
15 typedef Eigen::SparseMatrix<double> SpMat;
16 
17 using namespace std;
18 
19 static double sqr(double x) { return x * x; }
20 
21 namespace {
22 LOG_LOGGER _log = LOG_GET("jointcal.PhotometryFit");
23 }
24 
25 namespace lsst {
26 namespace jointcal {
27 
28 PhotometryFit::PhotometryFit(Associations &associations, PhotometryModel *photometryModel)
29  : _associations(associations),
30  _photometryModel(photometryModel),
31  _lastNTrip(0) {
32  // The various _npar... are initialized in assignIndices.
33  // Although there is no reason to adress them before one might be tempted by
34  // evaluating a Chi2 rightaway, .. which uses these counts, so:
35  assignIndices("");
36 }
37 
38 void PhotometryFit::LSDerivatives(TripletList &tripletList, Eigen::VectorXd &rhs) const {
39  auto ccdImageList = _associations.getCcdImageList();
40  for (auto const &im : ccdImageList) {
41  LSDerivatives(*im, tripletList, rhs);
42  }
43 }
44 
45 void PhotometryFit::LSDerivatives(const CcdImage &ccdImage, TripletList &tripletList, Eigen::VectorXd &rhs,
46  const MeasuredStarList *measuredStarList) const {
47  /***************************************************************************/
49  /***************************************************************************/
50  /* this routine works in two different ways: either providing the
51  Ccd, of providing the MeasuredStarList. In the latter case, the
52  Ccd should match the one(s) in the list. */
53  if (measuredStarList) assert(&(measuredStarList->front()->getCcdImage()) == &ccdImage);
54 
55  unsigned npar_max = 100; // anything large
56  vector<unsigned> indices(npar_max, -1);
57 
58  Eigen::VectorXd h(npar_max);
59  Eigen::VectorXd grad(npar_max);
60  // current position in the Jacobian
61  unsigned kTriplets = tripletList.getNextFreeIndex();
62  const MeasuredStarList &catalog = (measuredStarList) ? *measuredStarList : ccdImage.getCatalogForFit();
63 
64  for (auto const &i : catalog) {
65  const MeasuredStar &measuredStar = *i;
66  if (!measuredStar.isValid()) continue;
67  // tweak the measurement errors
68  double sigma = measuredStar.eflux;
69 #ifdef FUTURE
70  TweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
71 #endif
72  h.setZero(); // we cannot be sure that all entries will be overwritten.
73 
74  double pf = _photometryModel->photomFactor(ccdImage, measuredStar);
75  auto fs = measuredStar.getFittedStar();
76 
77  double res = measuredStar.getFlux() - pf * fs->getFlux();
78 
79  if (_fittingModel) {
80  _photometryModel->getIndicesAndDerivatives(measuredStar, ccdImage, indices, h);
81  for (unsigned k = 0; k < indices.size(); k++) {
82  unsigned l = indices[k];
83  tripletList.addTriplet(l, kTriplets, h[k] * fs->getFlux() / sigma);
84  rhs[l] += h[k] * res / sqr(sigma);
85  }
86  }
87  if (_fittingFluxes) {
88  unsigned index = fs->getIndexInMatrix();
89  tripletList.addTriplet(index, kTriplets, pf / sigma);
90  rhs[index] += res * pf / sqr(sigma);
91  }
92  kTriplets += 1; // each measurement contributes 1 column in the Jacobian
93  } // end loop on measurements
94  tripletList.setNextFreeIndex(kTriplets);
95 }
96 
97 // This is almost a selection of lines of LSDerivatives(CcdImage ...)
98 /* This routine is template because it is used
99 both with its first argument as "const CcdImageList &" and "CcdImageList &",
100 and I did not want to replicate it. The constness of the iterators is
101 automagically set by declaring them as "auto" */
102 
103 template <class ListType, class Accum>
104 void PhotometryFit::accumulateStat(ListType &listType, Accum &accum) const {
105  for (auto &im : listType) {
106  /**********************************************************************/
108  /**********************************************************************/
109  auto &ccdIMage = *im;
110  auto &catalog = ccdIMage.getCatalogForFit();
111 
112  for (auto const &measuredStar : catalog) {
113  if (!measuredStar->isValid()) continue;
114  // tweak the measurement errors
115  double sigma = measuredStar->eflux;
116 #ifdef FUTURE
117  TweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
118 #endif
119 
120  double pf = _photometryModel->photomFactor(ccdIMage, *measuredStar);
121  auto fs = measuredStar->getFittedStar();
122  double res = measuredStar->getFlux() - pf * fs->getFlux();
123  double chi2Val = sqr(res / sigma);
124  accum.addEntry(chi2Val, 1, measuredStar);
125  } // end loop on measurements
126  }
127 }
128 
131  Chi2 chi2;
132  accumulateStat(_associations.getCcdImageList(), chi2);
133  // so far, chi2.ndof contains the number of squares.
134  // So, subtract here the number of parameters.
135  chi2.ndof -= _nParTot;
136  return chi2;
137 }
138 
139 void PhotometryFit::outliersContributions(MeasuredStarList &outliers, TripletList &tripletList,
140  Eigen::VectorXd &grad) {
141  for (auto &outlier : outliers) {
142  MeasuredStarList tmp;
143  tmp.push_back(outlier);
144  const CcdImage &ccdImage = outlier->getCcdImage();
145  LSDerivatives(ccdImage, tripletList, grad, &tmp);
146  outlier->setValid(false);
147  auto fs = std::const_pointer_cast<FittedStar>(outlier->getFittedStar());
148  fs->getMeasurementCount()--;
149  }
150 }
151 
153 
157 struct Chi2Entry {
158  double chi2;
159  std::shared_ptr<MeasuredStar> measuredStar;
160 
161  Chi2Entry(double chi2, std::shared_ptr<MeasuredStar> star) : chi2(chi2), measuredStar(std::move(star)) {}
162  // for sort
163  bool operator<(const Chi2Entry &right) const { return (chi2 < right.chi2); }
164 };
165 
166 struct Chi2Vect : public vector<Chi2Entry> {
167  void addEntry(double Chi2Val, unsigned ndof, std::shared_ptr<MeasuredStar> measuredStar) {
168  push_back(Chi2Entry(Chi2Val, std::move(measuredStar)));
169  }
170 };
171 
173 
175 void PhotometryFit::getMeasuredStarIndices(const MeasuredStar &measuredStar,
176  std::vector<unsigned> &indices) const {
177  indices.clear();
178  if (_fittingModel) {
179  Eigen::VectorXd h(100);
180  _photometryModel->getIndicesAndDerivatives(measuredStar, measuredStar.getCcdImage(), indices, h);
181  }
182  if (_fittingFluxes) {
183  auto fs = measuredStar.getFittedStar();
184  unsigned fsIndex = fs->getIndexInMatrix();
185  indices.push_back(fsIndex);
186  }
187 }
188 
189 void PhotometryFit::findOutliers(double nSigCut, MeasuredStarList &outliers) const {
190  /* Aims at providing an outlier list for small-rank update
191  of the factorization. */
192 
193  // collect chi2 contributions of the measurements.
194  Chi2Vect chi2s;
195  // chi2s.reserve(_nMeasuredStars);
196  accumulateStat(_associations.ccdImageList, chi2s);
197  // do some stat
198  unsigned nval = chi2s.size();
199  if (nval == 0) return;
200  sort(chi2s.begin(), chi2s.end()); // increasing order. We rely on this later.
201  double median =
202  (nval & 1) ? chi2s[nval / 2].chi2 : 0.5 * (chi2s[nval / 2 - 1].chi2 + chi2s[nval / 2].chi2);
203  // some more stats. should go into the class if recycled anywhere else
204  double sum = 0;
205  double sum2 = 0;
206  for (auto const &i : chi2s) {
207  sum += i.chi2;
208  sum2 += sqr(i.chi2);
209  }
210  double average = sum / nval;
211  double sigma = sqrt(sum2 / nval - sqr(average));
212  LOGLS_INFO(_log,
213  "findOutliers chi2 stat: mean/median/sigma " << average << '/' << median << '/' << sigma);
214  double cut = average + nSigCut * sigma;
215  /* For each of the parameters, we will not remove more than 1
216  measurement that contributes to constraining it. Keep track
217  of the affected parameters using an integer vector. This is the
218  trick that Marc Betoule came up to for outlier removals in "star
219  flats" fits. */
220  Eigen::VectorXi affectedParams(_nParTot);
221  affectedParams.setZero();
222 
223  // start from the strongest outliers, i.e. at the end of the array.
224  for (auto i = chi2s.rbegin(); i != chi2s.rend(); ++i) {
225  if (i->chi2 < cut) break; // because the array is sorted.
226  vector<unsigned> indices;
227  getMeasuredStarIndices(*(i->measuredStar), indices);
228  bool drop_it = true;
229  /* find out if a stronger outlier constraining one of the parameters
230  this one contrains was already discarded. If yes, we keep this one */
231  for (auto const &i : indices)
232  if (affectedParams(i) != 0) drop_it = false;
233 
234  if (drop_it) {
235  for (auto const &i : indices) affectedParams(i)++;
236  outliers.push_back(i->measuredStar);
237  }
238  } // end loop on measurements
239  LOGLS_INFO(_log, "findMeasOutliers: found " << outliers.size() << " outliers");
240 }
241 
242 void PhotometryFit::assignIndices(const std::string &whatToFit) {
243  _whatToFit = whatToFit;
244  LOGLS_INFO(_log, "assignIndices: now fitting: " << whatToFit);
245  _fittingModel = (_whatToFit.find("Model") != string::npos);
246  _fittingFluxes = (_whatToFit.find("Fluxes") != string::npos);
247  // When entering here, we assume that whatToFit has already been interpreted.
248 
249  _nParModel = (_fittingModel) ? _photometryModel->assignIndices(whatToFit, 0) : 0;
250  unsigned ipar = _nParModel;
251 
252  if (_fittingFluxes) {
253  FittedStarList &fsl = _associations.fittedStarList;
254  for (auto &i : fsl) {
255  FittedStar &fs = *i;
256  // the parameter layout here is used also
257  // - when filling the derivatives
258  // - when updating (offsetParams())
259  // - in GetMeasuredStarIndices
260  fs.setIndexInMatrix(ipar);
261  ipar += 1;
262  }
263  }
264  _nParTot = ipar;
265 }
266 
267 void PhotometryFit::offsetParams(const Eigen::VectorXd &delta) {
268  if (delta.size() != _nParTot)
269  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
270  "PhotometryFit::offsetParams : the provided vector length is not compatible with "
271  "the current whatToFit setting");
272  if (_fittingModel) _photometryModel->offsetParams(delta);
273 
274  if (_fittingFluxes) {
275  FittedStarList &fsl = _associations.fittedStarList;
276  for (auto &i : fsl) {
277  FittedStar &fs = *i;
278  // the parameter layout here is used also
279  // - when filling the derivatives
280  // - when assigning indices (assignIndices())
281  unsigned index = fs.getIndexInMatrix();
282  fs.getFlux() += delta(index);
283  }
284  }
285 }
286 
287 bool PhotometryFit::minimize(const std::string &whatToFit) {
288  assignIndices(whatToFit);
289 
290  // TODO : write a guesser for the number of triplets
291  unsigned nTrip = (_lastNTrip) ? _lastNTrip : 1e6;
292  TripletList tripletList(nTrip);
293  Eigen::VectorXd grad(_nParTot);
294  grad.setZero();
295 
296  // Fill the triplets
297  LSDerivatives(tripletList, grad);
298  _lastNTrip = tripletList.size();
299 
300  SpMat hessian;
301  {
302  SpMat jacobian(_nParTot, tripletList.getNextFreeIndex());
303  jacobian.setFromTriplets(tripletList.begin(), tripletList.end());
304  // release memory shrink_to_fit is C++11
305  tripletList.clear(); // tripletList.shrink_to_fit();
306  hessian = jacobian * jacobian.transpose();
307  } // release the Jacobian
308 
309  LOGLS_DEBUG(_log, "Starting factorization, hessian: dim="
310  << hessian.rows() << " nnz=" << hessian.nonZeros()
311  << " filling-frac = " << hessian.nonZeros() / sqr(hessian.rows()));
312 
313  Eigen::SimplicialLDLT<SpMat> chol(hessian);
314  if (chol.info() != Eigen::Success) {
315  LOGLS_ERROR(_log, "minimize: factorization failed ");
316  return false;
317  }
318 
319  Eigen::VectorXd delta = chol.solve(grad);
320 
321  offsetParams(delta);
322  return true;
323 }
324 
325 void PhotometryFit::makeResTuple(const std::string &tupleName) const {
326  std::ofstream tuple(tupleName.c_str());
327  /* If we think the some coordinate on the focal plane is relevant in
328  the ntuple, because thmodel relies on it, then we have to add
329  some function to the model that returns this relevant
330  coordinate. */
331  tuple << "#xccd: coordinate in CCD" << endl
332  << "#yccd: " << endl
333  << "#mag: rough mag" << endl
334  << "#flux : measured flux" << endl
335  << "#eflux : measured flux erro" << endl
336  << "#fflux : fitted flux" << endl
337  << "#phot_factor:" << endl
338  << "#jd: Julian date of the measurement" << endl
339  << "#color : " << endl
340  << "#fsindex: some unique index of the object" << endl
341  << "#ra: pos of fitted star" << endl
342  << "#dec: pos of fitted star" << endl
343  << "#chi2: contribution to Chi2 (1 dof)" << endl
344  << "#nm: number of measurements of this FittedStar" << endl
345  << "#chip: chip number" << endl
346  << "#visit: visit id" << endl
347  << "#end" << endl;
348  const CcdImageList &ccdImageList = _associations.getCcdImageList();
349  for (auto const &i : ccdImageList) {
350  const CcdImage &im = *i;
351  const MeasuredStarList &cat = im.getCatalogForFit();
352  for (auto const &is : cat) {
353  const MeasuredStar &ms = *is;
354  if (!ms.isValid()) continue;
355  double sigma = ms.eflux;
356 #ifdef FUTURE
357  tweakPhotomMeasurementErrors(inPos, ms, _fluxError);
358 #endif
359  double pf = _photometryModel->photomFactor(im, ms);
360  double jd = im.getMjd();
361  auto fs = ms.getFittedStar();
362  double res = ms.getFlux() - pf * fs->getFlux();
363  double chi2Val = sqr(res / sigma);
364  tuple << ms.x << ' ' << ms.y << ' ' << fs->getMag() << ' ' << ms.getFlux() << ' ' << ms.eflux
365  << ' ' << fs->getFlux() << ' ' << pf << ' ' << jd << ' ' << fs->color << ' '
366  << fs->getIndexInMatrix() << ' ' << fs->x << ' ' << fs->y << ' ' << chi2Val << ' '
367  << fs->getMeasurementCount() << ' ' << im.getCcdId() << ' ' << im.getVisit() << endl;
368  } // loop on measurements in image
369  } // loop on images
370 }
371 } // namespace jointcal
372 } // namespace lsst
Chi2 computeChi2() const
Returns a chi2 for the current state.
double getMjd() const
Julian Date.
Definition: CcdImage.h:136
a class to accumulate chi2 contributions together with pointers to the contributors.
void LSDerivatives(TripletList &tripletList, Eigen::VectorXd &rhs) const
Derivatives of the Chi2.
bool minimize(const std::string &whatToFit)
Does a 1 step minimization, assuming a linear model.
virtual double photomFactor(const CcdImage &ccdImage, const Point &where) const =0
Return the &quot;photometric factor&quot; at a given location on a ccdImage.
std::shared_ptr< const FittedStar > getFittedStar() const
Definition: MeasuredStar.h:56
void addEntry(double Chi2Val, unsigned ndof, std::shared_ptr< MeasuredStar > measuredStar)
double chi2
The class that implements the relations between MeasuredStar and FittedStar.
Definition: Associations.h:28
unsigned ndof
Definition: Chi2.h:15
A list of MeasuredStar. They are usually filled in Associations::AddImage.
Definition: MeasuredStar.h:71
void assignIndices(const std::string &whatToFit)
Set parameter groups fixed or variable and assign indices to each parameter in the big matrix (which ...
virtual void offsetParams(const Eigen::VectorXd &delta)=0
Offset the parameters by the provided amounts.
double x
coordinate
Definition: Point.h:18
Simple structure to accumulate Chi2 and Ndof.
Definition: Chi2.h:13
const CcdImageList & getCcdImageList() const
Definition: Associations.h:109
bool operator<(const Chi2Entry &right) const
double getFlux() const
Definition: BaseStar.h:71
void setIndexInMatrix(const unsigned &index)
index is a value that a fit can set and reread....
Definition: FittedStar.h:114
int getIndexInMatrix() const
Definition: FittedStar.h:117
void addTriplet(const unsigned i, const unsigned j, double val)
Definition: Tripletlist.h:23
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition: FittedStar.h:134
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
virtual void getIndicesAndDerivatives(const MeasuredStar &measuredStar, const CcdImage &ccdImage, std::vector< unsigned > &indices, Eigen::VectorXd &D)=0
number of parameters to be read in indices.size()
unsigned getNextFreeIndex() const
Definition: Tripletlist.h:25
Interface class for PhotometryFit.
Chi2Entry(double chi2, std::shared_ptr< MeasuredStar > star)
int getMeasurementCount() const
Definition: FittedStar.h:98
std::shared_ptr< MeasuredStar > measuredStar
virtual unsigned assignIndices(const std::string &whatToFit, unsigned firstIndex)=0
Assign indices to parameters involved in mappings, starting at firstIndex.
std::list< std::shared_ptr< CcdImage > > CcdImageList
Definition: CcdImage.h:22
void makeResTuple(const std::string &tupleName) const
Produces an ntuple.
Handler of an actual image from a single CCD.
Definition: CcdImage.h:31
Eigen::SparseMatrix< double > SpMat
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.
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
int getCcdId() const
returns ccd ID
Definition: CcdImage.h:127
FittedStarList fittedStarList
Definition: Associations.h:32