lsst.jointcal  21.0.0-19-g8edf490+4be6a9e008
FitterBase.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 /*
3  * This file is part of jointcal.
4  *
5  * Developed for the LSST Data Management System.
6  * This product includes software developed by the LSST Project
7  * (https://www.lsst.org).
8  * See the COPYRIGHT file at the top-level directory of this distribution
9  * for details of code ownership.
10  *
11  * This program is free software: you can redistribute it and/or modify
12  * it under the terms of the GNU General Public License as published by
13  * the Free Software Foundation, either version 3 of the License, or
14  * (at your option) any later version.
15  *
16  * This program is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19  * GNU General Public License for more details.
20  *
21  * You should have received a copy of the GNU General Public License
22  * along with this program. If not, see <https://www.gnu.org/licenses/>.
23  */
24 
25 #include <vector>
26 #include "Eigen/Core"
27 
28 #include <boost/math/tools/minima.hpp>
29 
30 #include "lsst/log/Log.h"
31 
32 #include "lsst/jointcal/Chi2.h"
33 #include "lsst/jointcal/CcdImage.h"
38 
39 namespace lsst {
40 namespace jointcal {
41 
43  Chi2Statistic chi2;
44  accumulateStatImageList(_associations->getCcdImageList(), chi2);
46  // chi2.ndof contains the number of squares.
47  // So subtract the number of parameters.
48  chi2.ndof -= _nParTot;
49  return chi2;
50 }
51 
53  FittedStarList &fsOutliers, double &cut) const {
54  // collect chi2 contributions
55  Chi2List chi2List;
56  chi2List.reserve(_nMeasuredStars + _associations->refStarList.size());
57  // contributions from measurement terms:
58  accumulateStatImageList(_associations->ccdImageList, chi2List);
59  // and from reference terms
60  accumulateStatRefStars(chi2List);
61 
62  // compute some statistics
63  size_t nval = chi2List.size();
64  if (nval == 0) return 0;
65  sort(chi2List.begin(), chi2List.end());
66  double median = (nval & 1) ? chi2List[nval / 2].chi2
67  : 0.5 * (chi2List[nval / 2 - 1].chi2 + chi2List[nval / 2].chi2);
68  auto averageAndSigma = chi2List.computeAverageAndSigma();
69  LOGLS_DEBUG(_log, "findOutliers chi2 stat: mean/median/sigma " << averageAndSigma.first << '/' << median
70  << '/' << averageAndSigma.second);
71  cut = averageAndSigma.first + nSigmaCut * averageAndSigma.second;
72  /* For each of the parameters, we will not remove more than 1
73  measurement that contributes to constraining it. Keep track using
74  of what we are touching using an integer vector. This is the
75  trick that Marc Betoule came up to for outlier removals in "star
76  flats" fits. */
77  Eigen::VectorXi affectedParams(_nParTot);
78  affectedParams.setZero();
79 
80  std::size_t nOutliers = 0; // returned to the caller
81  // start from the strongest outliers.
82  for (auto chi2 = chi2List.rbegin(); chi2 != chi2List.rend(); ++chi2) {
83  if (chi2->chi2 < cut) break; // because the array is sorted.
84  IndexVector indices;
85  /* now, we want to get the indices of the parameters this chi2
86  term depends on. We have to figure out which kind of term it
87  is; we use for that the type of the star attached to the Chi2Star. */
88  auto measuredStar = std::dynamic_pointer_cast<MeasuredStar>(chi2->star);
89  std::shared_ptr<FittedStar> fittedStar; // To add to fsOutliers if it is a reference outlier.
90  if (measuredStar == nullptr) {
91  // it is a reference outlier
92  fittedStar = std::dynamic_pointer_cast<FittedStar>(chi2->star);
93  if (fittedStar->getMeasurementCount() == 0) {
94  LOGLS_WARN(_log, "FittedStar with no measuredStars found as an outlier: " << *fittedStar);
95  continue;
96  }
97  // NOTE: Stars contribute twice to astrometry (x,y), but once to photometry (flux),
98  // NOTE: but we only need to mark one index here because both will be removed with that star.
99  indices.push_back(fittedStar->getIndexInMatrix());
100  LOGLS_TRACE(_log, "Removing refStar " << *(fittedStar->getRefStar()) << " chi2: " << chi2->chi2);
101  /* One might think it would be useful to account for PM
102  parameters here, but it is just useless */
103  } else {
104  // it is a measurement outlier
105  auto tempFittedStar = measuredStar->getFittedStar();
106  if (tempFittedStar->getMeasurementCount() == 1 && tempFittedStar->getRefStar() == nullptr) {
107  LOGLS_WARN(_log, "FittedStar with 1 measuredStar and no refStar found as an outlier: "
108  << *tempFittedStar);
109  continue;
110  }
111  getIndicesOfMeasuredStar(*measuredStar, indices);
112  LOGLS_TRACE(_log, "Removing measStar " << *measuredStar << " chi2: " << chi2->chi2);
113  }
114 
115  /* Find out if we already discarded a stronger outlier
116  constraining some parameter this one constrains as well. If
117  yes, we keep this one, because this stronger outlier could be
118  causing the large chi2 we have in hand. */
119  bool drop_it = true;
120  for (auto const &i : indices) {
121  if (affectedParams(i) != 0) {
122  drop_it = false;
123  }
124  }
125 
126  if (drop_it) // store the outlier in one of the lists:
127  {
128  if (measuredStar == nullptr) {
129  // reference term
130  fsOutliers.push_back(fittedStar);
131  } else {
132  // measurement term
133  msOutliers.push_back(measuredStar);
134  }
135  // mark the parameters as directly changed when we discard this chi2 term.
136  for (auto const &i : indices) {
137  affectedParams(i)++;
138  }
139  nOutliers++;
140  }
141  } // end loop on measurements/references
142  LOGLS_INFO(_log, "findOutliers: found " << msOutliers.size() << " meas outliers and " << fsOutliers.size()
143  << " ref outliers ");
144 
145  return nOutliers;
146 }
147 
148 namespace {
150 SparseMatrixD createHessian(std::size_t nParTot, TripletList const &tripletList) {
151  SparseMatrixD jacobian(nParTot, tripletList.getNextFreeIndex());
152  jacobian.setFromTriplets(tripletList.begin(), tripletList.end());
153  return jacobian * jacobian.transpose();
154 }
155 
157 void dumpMatrixAndGradient(SparseMatrixD const &matrix, Eigen::VectorXd const &grad,
158  std::string const &dumpFile, LOG_LOGGER _log) {
159  std::string ext = ".txt";
160  Eigen::MatrixXd matrixDense(matrix);
161  std::string dumpMatrixPath = dumpFile + "-mat" + ext;
162  std::ofstream matfile(dumpMatrixPath);
163  matfile << matrixDense << std::endl;
164  std::string dumpGradPath = dumpFile + "-grad" + ext;
165  std::ofstream gradfile(dumpGradPath);
166  gradfile << grad << std::endl;
167  LOGLS_INFO(_log, "Dumped Hessian, gradient to: '" << dumpMatrixPath << "', '" << dumpGradPath << "'");
168 }
169 } // namespace
170 
171 MinimizeResult FitterBase::minimize(std::string const &whatToFit, double nSigmaCut,
172  double sigmaRelativeTolerance, bool doRankUpdate, bool const doLineSearch,
173  std::string const &dumpMatrixFile) {
174  assignIndices(whatToFit);
175 
177 
178  // TODO : write a guesser for the number of triplets
179  std::size_t nTrip = (_lastNTrip) ? _lastNTrip : 1e6;
180  TripletList tripletList(nTrip);
181  Eigen::VectorXd grad(_nParTot);
182  grad.setZero();
183  double scale = 1.0;
184 
185  // Fill the triplets
186  leastSquareDerivatives(tripletList, grad);
187  _lastNTrip = tripletList.size();
188 
189  LOGLS_DEBUG(_log, "End of triplet filling, ntrip = " << tripletList.size());
190 
191  SparseMatrixD hessian = createHessian(_nParTot, tripletList);
192  tripletList.clear(); // we don't need it any more after we have the hessian.
193 
194  LOGLS_DEBUG(_log, "Starting factorization, hessian: dim="
195  << hessian.rows() << " non-zeros=" << hessian.nonZeros()
196  << " filling-frac = " << hessian.nonZeros() / std::pow(hessian.rows(), 2));
197 
198  if (dumpMatrixFile != "") {
199  if (hessian.rows() * hessian.cols() > 2e8) {
200  LOGLS_WARN(_log, "Hessian matrix is too big to dump to file, with rows, columns: "
201  << hessian.rows() << ", " << hessian.cols());
202  } else {
203  dumpMatrixAndGradient(hessian, grad, dumpMatrixFile, _log);
204  }
205  }
206 
208  if (chol.info() != Eigen::Success) {
209  LOGLS_ERROR(_log, "minimize: factorization failed ");
210  return MinimizeResult::Failed;
211  }
212 
213  std::size_t totalMeasOutliers = 0;
214  std::size_t totalRefOutliers = 0;
215  double oldChi2 = computeChi2().chi2;
216  double oldSigmaCut = 0;
217  double sigmaCut;
218 
219  while (true) {
220  Eigen::VectorXd delta = chol.solve(grad);
221  if (doLineSearch) {
222  scale = _lineSearch(delta);
223  }
224  offsetParams(scale * delta);
225  Chi2Statistic currentChi2(computeChi2());
226  LOGLS_DEBUG(_log, currentChi2);
227  if (!isfinite(currentChi2.chi2)) {
228  LOGL_ERROR(_log, "chi2 is not finite. Aborting outlier rejection.");
229  returnCode = MinimizeResult::NonFinite;
230  break;
231  }
232  if (currentChi2.chi2 > oldChi2 && totalMeasOutliers + totalRefOutliers != 0) {
233  LOGL_WARN(_log, "chi2 went up, skipping outlier rejection loop");
234  returnCode = MinimizeResult::Chi2Increased;
235  break;
236  }
237  oldChi2 = currentChi2.chi2;
238 
239  if (nSigmaCut == 0) break; // no rejection step to perform
240  MeasuredStarList msOutliers;
241  FittedStarList fsOutliers;
242  // keep nOutliers so we don't have to sum msOutliers.size()+fsOutliers.size() twice below.
243  std::size_t nOutliers = findOutliers(nSigmaCut, msOutliers, fsOutliers, sigmaCut);
244  double relChange = 1 - sigmaCut / oldSigmaCut;
245  LOGLS_DEBUG(_log, "findOutliers chi2 cut level: " << sigmaCut << ", relative change: " << relChange);
246  // If sigmaRelativeTolerance is set and at least one iteration has been done, break loop when the
247  // fractional change in sigmaCut levels is less than the sigmaRelativeTolerance parameter.
248  if ((sigmaRelativeTolerance > 0) && (oldSigmaCut > 0 && relChange < sigmaRelativeTolerance)){
249  LOGLS_INFO(_log, "Iterations stopped with chi2 cut at " << sigmaCut << " and relative change of "
250  << relChange);
251  break;
252  }
253  totalMeasOutliers += msOutliers.size();
254  totalRefOutliers += fsOutliers.size();
255  oldSigmaCut = sigmaCut;
256  if (nOutliers == 0) break;
257  TripletList outlierTriplets(nOutliers);
258  grad.setZero(); // recycle the gradient
259  // compute the contributions of outliers to derivatives
260  outliersContributions(msOutliers, fsOutliers, outlierTriplets, grad);
261  // Remove significant outliers
262  removeMeasOutliers(msOutliers);
263  removeRefOutliers(fsOutliers);
264  if (doRankUpdate) {
265  // convert triplet list to eigen internal format
266  SparseMatrixD H(_nParTot, outlierTriplets.getNextFreeIndex());
267  H.setFromTriplets(outlierTriplets.begin(), outlierTriplets.end());
268  chol.update(H, false /* means downdate */);
269  // The contribution of outliers to the gradient is the opposite
270  // of the contribution of all other terms, because they add up to 0
271  grad *= -1;
272  } else {
273  // don't reuse tripletList because we want a new nextFreeIndex.
274  TripletList nextTripletList(_lastNTrip);
275  grad.setZero();
276  // Rebuild the matrix and gradient
277  leastSquareDerivatives(nextTripletList, grad);
278  _lastNTrip = nextTripletList.size();
279  LOGLS_DEBUG(_log, "Triplets recomputed, ntrip = " << nextTripletList.size());
280 
281  hessian = createHessian(_nParTot, nextTripletList);
282  nextTripletList.clear(); // we don't need it any more after we have the hessian.
283 
285  "Restarting factorization, hessian: dim="
286  << hessian.rows() << " non-zeros=" << hessian.nonZeros()
287  << " filling-frac = " << hessian.nonZeros() / std::pow(hessian.rows(), 2));
288  chol.compute(hessian);
289  if (chol.info() != Eigen::Success) {
290  LOGLS_ERROR(_log, "minimize: factorization failed ");
291  return MinimizeResult::Failed;
292  }
293  }
294  }
295 
296  // only print the outlier summary if outlier rejection was turned on.
297  if (nSigmaCut != 0) {
298  LOGLS_INFO(_log, "Number of outliers (Measured + Reference = Total): "
299  << totalMeasOutliers << " + " << totalRefOutliers << " = "
300  << totalMeasOutliers + totalRefOutliers);
301  }
302  return returnCode;
303 }
304 
306  TripletList &tripletList, Eigen::VectorXd &grad) {
307  for (auto &outlier : msOutliers) {
308  MeasuredStarList tmp;
309  tmp.push_back(outlier);
310  const CcdImage &ccdImage = outlier->getCcdImage();
311  leastSquareDerivativesMeasurement(ccdImage, tripletList, grad, &tmp);
312  }
313  leastSquareDerivativesReference(fsOutliers, tripletList, grad);
314 }
315 
317  for (auto &measuredStar : outliers) {
318  auto fittedStar = measuredStar->getFittedStar();
319  measuredStar->setValid(false);
320  fittedStar->getMeasurementCount()--; // could be put in setValid
321  }
322 }
323 
325  for (auto &fittedStar : outliers) {
326  fittedStar->setRefStar(nullptr);
327  }
328 }
329 
330 void FitterBase::leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const {
331  auto ccdImageList = _associations->getCcdImageList();
332  for (auto const &ccdImage : ccdImageList) {
333  leastSquareDerivativesMeasurement(*ccdImage, tripletList, grad);
334  }
335  leastSquareDerivativesReference(_associations->fittedStarList, tripletList, grad);
336 }
337 
338 void FitterBase::saveChi2Contributions(std::string const &baseName) const {
339  std::string replaceStr = "{type}";
340  auto pos = baseName.find(replaceStr);
341  std::string measFilename(baseName);
342  measFilename.replace(pos, replaceStr.size(), "-meas.csv");
343  std::string refFilename(baseName);
344  refFilename.replace(pos, replaceStr.size(), "-ref.csv");
345  saveChi2MeasContributions(measFilename);
346  saveChi2RefContributions(refFilename);
347 }
348 
349 double FitterBase::_lineSearch(Eigen::VectorXd const &delta) {
350  auto func = [this, &delta](double scale) {
351  auto offset = scale * delta;
352  offsetParams(offset);
353  auto chi2 = computeChi2();
354  // reset the system to where it was before offsetting.
355  offsetParams(-offset);
356  return chi2.chi2;
357  };
358  // The maximum theoretical precision is half the number of bits in the mantissa (see boost docs).
359  auto bits = std::numeric_limits<double>::digits / 2;
360  auto result = boost::math::tools::brent_find_minima(func, -1.0, 2.0, bits);
361  LOGLS_DEBUG(_log, "Line search scale factor: " << result.first);
362  return result.first;
363 }
364 
365 } // namespace jointcal
366 } // namespace lsst
Eigen::SparseMatrix< double, 0, Eigen::Index > SparseMatrixD
Definition: Eigenstuff.h:35
#define LOGLS_WARN(logger, message)
#define LOGL_WARN(logger, message...)
#define LOGLS_INFO(logger, message)
#define LOGLS_ERROR(logger, message)
#define LOGL_ERROR(logger, message...)
#define LOGLS_DEBUG(logger, message)
#define LOGLS_TRACE(logger, message)
T begin(T... args)
void update(SparseMatrixD const &H, bool UpOrDown)
Definition: Eigenstuff.h:68
Handler of an actual image from a single CCD.
Definition: CcdImage.h:64
Structure to accumulate the chi2 contributions per each star (to help find outliers).
Definition: Chi2.h:100
std::pair< double, double > computeAverageAndSigma()
Compute the average and std-deviation of these chisq values.
Definition: Chi2.cc:33
Simple structure to accumulate chi2 and ndof.
Definition: Chi2.h:52
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition: FittedStar.h:123
void leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const
Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting.
Definition: FitterBase.cc:330
void removeRefOutliers(FittedStarList &outliers)
Remove refStar outliers from the fit. No Refit done.
Definition: FitterBase.cc:324
virtual void getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, IndexVector &indices) const =0
Set the indices of a measured star from the full matrix, for outlier removal.
Chi2Statistic computeChi2() const
Returns the chi2 for the current state.
Definition: FitterBase.cc:42
virtual void saveChi2MeasContributions(std::string const &filename) const =0
Save a CSV file containing residuals of measurement terms.
virtual void leastSquareDerivativesReference(FittedStarList const &fittedStarList, TripletList &tripletList, Eigen::VectorXd &grad) const =0
Compute the derivatives of the reference terms.
virtual void saveChi2Contributions(std::string const &baseName) const
Save the full chi2 term per star that was used in the minimization, for debugging.
Definition: FitterBase.cc:338
Eigen::Index _nMeasuredStars
Definition: FitterBase.h:161
virtual void assignIndices(std::string const &whatToFit)=0
Set parameters to fit and assign indices in the big matrix.
virtual void offsetParams(Eigen::VectorXd const &delta)=0
Offset the parameters by the requested quantities.
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:156
virtual void accumulateStatRefStars(Chi2Accumulator &accum) const =0
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars.
virtual void saveChi2RefContributions(std::string const &filename) const =0
Save a CSV file containing residuals of reference terms.
virtual void leastSquareDerivativesMeasurement(CcdImage const &ccdImage, TripletList &tripletList, Eigen::VectorXd &grad, MeasuredStarList const *measuredStarList=nullptr) const =0
Compute the derivatives of the measured stars and model for one CcdImage.
void outliersContributions(MeasuredStarList &msOutliers, FittedStarList &fsOutliers, TripletList &tripletList, Eigen::VectorXd &grad)
Contributions to derivatives from (presumably) outlier terms.
Definition: FitterBase.cc:305
std::size_t findOutliers(double nSigmaCut, MeasuredStarList &msOutliers, FittedStarList &fsOutliers, double &cut) const
Find Measurements and references contributing more than a cut, computed as.
Definition: FitterBase.cc:52
virtual void accumulateStatImageList(CcdImageList const &ccdImageList, Chi2Accumulator &accum) const =0
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for measurements.
void removeMeasOutliers(MeasuredStarList &outliers)
Remove measuredStar outliers from the fit. No Refit done.
Definition: FitterBase.cc:316
MinimizeResult minimize(std::string const &whatToFit, double const nSigmaCut=0, double sigmaRelativeTolerance=0, bool const doRankUpdate=true, bool const doLineSearch=false, std::string const &dumpMatrixFile="")
Does a 1 step minimization, assuming a linear model.
Definition: FitterBase.cc:171
A list of MeasuredStar. They are usually filled in Associations::createCcdImage.
Definition: MeasuredStar.h:146
Eigen::Index getNextFreeIndex() const
Definition: Tripletlist.h:47
T clear(T... args)
T end(T... args)
T endl(T... args)
T find(T... args)
T isfinite(T... args)
def scale(algorithm, min, max=None, frame=None)
MinimizeResult
Return value of minimize()
Definition: FitterBase.h:40
Class for a simple mapping implementing a generic AstrometryTransform.
T pow(T... args)
T push_back(T... args)
T rbegin(T... args)
T rend(T... args)
T replace(T... args)
T reserve(T... args)
T size(T... args)
T sort(T... args)