lsst.jointcal  15.0-17-g076ea75+10
FitterBase.cc
Go to the documentation of this file.
1 #include <vector>
2 #include "Eigen/Core"
3 
4 #include "lsst/log/Log.h"
5 
6 #include "lsst/jointcal/Chi2.h"
12 
13 namespace lsst {
14 namespace jointcal {
15 
17  Chi2Statistic chi2;
18  accumulateStatImageList(_associations->getCcdImageList(), chi2);
20  // chi2.ndof contains the number of squares.
21  // So subtract the number of parameters.
22  chi2.ndof -= _nParTot;
23  return chi2;
24 }
25 
26 unsigned FitterBase::findOutliers(double nSigmaCut, MeasuredStarList &msOutliers,
27  FittedStarList &fsOutliers) const {
28  // collect chi2 contributions
29  Chi2List chi2List;
30  chi2List.reserve(_nMeasuredStars + _associations->refStarList.size());
31  // contributions from measurement terms:
32  accumulateStatImageList(_associations->ccdImageList, chi2List);
33  // and from reference terms
34  accumulateStatRefStars(chi2List);
35 
36  // compute some statistics
37  size_t nval = chi2List.size();
38  if (nval == 0) return 0;
39  sort(chi2List.begin(), chi2List.end());
40  double median = (nval & 1) ? chi2List[nval / 2].chi2
41  : 0.5 * (chi2List[nval / 2 - 1].chi2 + chi2List[nval / 2].chi2);
42  auto averageAndSigma = chi2List.computeAverageAndSigma();
43  LOGLS_DEBUG(_log, "RemoveOutliers chi2 stat: mean/median/sigma " << averageAndSigma.first << '/' << median
44  << '/' << averageAndSigma.second);
45  double cut = averageAndSigma.first + nSigmaCut * averageAndSigma.second;
46  /* For each of the parameters, we will not remove more than 1
47  measurement that contributes to constraining it. Keep track using
48  of what we are touching using an integer vector. This is the
49  trick that Marc Betoule came up to for outlier removals in "star
50  flats" fits. */
51  Eigen::VectorXi affectedParams(_nParTot);
52  affectedParams.setZero();
53 
54  unsigned nOutliers = 0; // returned to the caller
55  // start from the strongest outliers.
56  for (auto chi2 = chi2List.rbegin(); chi2 != chi2List.rend(); ++chi2) {
57  if (chi2->chi2 < cut) break; // because the array is sorted.
58  std::vector<unsigned> indices;
59  /* now, we want to get the indices of the parameters this chi2
60  term depends on. We have to figure out which kind of term it
61  is; we use for that the type of the star attached to the Chi2Star. */
62  auto measuredStar = std::dynamic_pointer_cast<MeasuredStar>(chi2->star);
63  std::shared_ptr<FittedStar> fittedStar; // To add to fsOutliers if it is a reference outlier.
64  if (measuredStar == nullptr) {
65  // it is a reference outlier
66  fittedStar = std::dynamic_pointer_cast<FittedStar>(chi2->star);
67  if (fittedStar->getMeasurementCount() == 0) {
68  LOGLS_WARN(_log, "FittedStar with no measuredStars found as an outlier: " << *fittedStar);
69  continue;
70  }
71  // NOTE: Stars contribute twice to astrometry (x,y), but once to photometry (flux),
72  // NOTE: but we only need to mark one index here because both will be removed with that star.
73  indices.push_back(fittedStar->getIndexInMatrix());
74  /* One might think it would be useful to account for PM
75  parameters here, but it is just useless */
76  } else {
77  // it is a measurement outlier
78  auto tempFittedStar = measuredStar->getFittedStar();
79  if (tempFittedStar->getMeasurementCount() == 1 && tempFittedStar->getRefStar() == nullptr) {
80  LOGLS_WARN(_log, "FittedStar with 1 measuredStar and no refStar found as an outlier: "
81  << *tempFittedStar);
82  continue;
83  }
84  getIndicesOfMeasuredStar(*measuredStar, indices);
85  }
86 
87  /* Find out if we already discarded a stronger outlier
88  constraining some parameter this one constrains as well. If
89  yes, we keep this one, because this stronger outlier could be
90  causing the large chi2 we have in hand. */
91  bool drop_it = true;
92  for (auto const &i : indices) {
93  if (affectedParams(i) != 0) {
94  drop_it = false;
95  }
96  }
97 
98  if (drop_it) // store the outlier in one of the lists:
99  {
100  if (measuredStar == nullptr) {
101  // reference term
102  fsOutliers.push_back(fittedStar);
103  } else {
104  // measurement term
105  msOutliers.push_back(measuredStar);
106  }
107  // mark the parameters as directly changed when we discard this chi2 term.
108  for (auto const &i : indices) {
109  affectedParams(i)++;
110  }
111  nOutliers++;
112  }
113  } // end loop on measurements/references
114  LOGLS_INFO(_log, "findOutliers: found " << msOutliers.size() << " meas outliers and " << fsOutliers.size()
115  << " ref outliers ");
116 
117  return nOutliers;
118 }
119 
120 namespace {
122 SparseMatrixD createHessian(int nParTot, TripletList const &tripletList) {
123  SparseMatrixD jacobian(nParTot, tripletList.getNextFreeIndex());
124  jacobian.setFromTriplets(tripletList.begin(), tripletList.end());
125  return jacobian * jacobian.transpose();
126 }
127 
129 void dumpMatrixAndGradient(SparseMatrixD const &matrix, Eigen::VectorXd const &grad,
130  std::string const &dumpFile, LOG_LOGGER _log) {
131  std::string ext = ".txt";
132  Eigen::MatrixXd matrixDense(matrix);
133  std::string dumpMatrixPath = dumpFile + "-mat" + ext;
134  std::ofstream matfile(dumpMatrixPath);
135  matfile << matrixDense << std::endl;
136  std::string dumpGradPath = dumpFile + "-grad" + ext;
137  std::ofstream gradfile(dumpGradPath);
138  gradfile << grad << std::endl;
139  LOGLS_INFO(_log, "Dumped Hessian, gradient to: '" << dumpMatrixPath << "', '" << dumpGradPath << "'");
140 }
141 } // namespace
142 
143 MinimizeResult FitterBase::minimize(std::string const &whatToFit, double nSigmaCut, bool doRankUpdate,
144  std::string const &dumpMatrixFile) {
145  assignIndices(whatToFit);
146 
148 
149  // TODO : write a guesser for the number of triplets
150  unsigned nTrip = (_lastNTrip) ? _lastNTrip : 1e6;
151  TripletList tripletList(nTrip);
152  Eigen::VectorXd grad(_nParTot);
153  grad.setZero();
154 
155  // Fill the triplets
156  leastSquareDerivatives(tripletList, grad);
157  _lastNTrip = tripletList.size();
158 
159  LOGLS_DEBUG(_log, "End of triplet filling, ntrip = " << tripletList.size());
160 
161  SparseMatrixD hessian = createHessian(_nParTot, tripletList);
162  tripletList.clear(); // we don't need it any more after we have the hessian.
163 
164  LOGLS_DEBUG(_log, "Starting factorization, hessian: dim="
165  << hessian.rows() << " non-zeros=" << hessian.nonZeros()
166  << " filling-frac = " << hessian.nonZeros() / std::pow(hessian.rows(), 2));
167 
168  if (dumpMatrixFile != "") {
169  if (hessian.rows() * hessian.cols() > 2e8) {
170  LOGLS_WARN(_log, "Hessian matrix is too big to dump to file, with rows, columns: "
171  << hessian.rows() << ", " << hessian.cols());
172  } else {
173  dumpMatrixAndGradient(hessian, grad, dumpMatrixFile, _log);
174  }
175  }
176 
178  if (chol.info() != Eigen::Success) {
179  LOGLS_ERROR(_log, "minimize: factorization failed ");
180  return MinimizeResult::Failed;
181  }
182 
183  unsigned totalMeasOutliers = 0;
184  unsigned totalRefOutliers = 0;
185  double oldChi2 = computeChi2().chi2;
186 
187  while (true) {
188  Eigen::VectorXd delta = chol.solve(grad);
189  offsetParams(delta);
190  Chi2Statistic currentChi2(computeChi2());
191  LOGLS_DEBUG(_log, currentChi2);
192  if (currentChi2.chi2 > oldChi2 && totalMeasOutliers + totalRefOutliers != 0) {
193  LOGL_WARN(_log, "chi2 went up, skipping outlier rejection loop");
194  returnCode = MinimizeResult::Chi2Increased;
195  break;
196  }
197  oldChi2 = currentChi2.chi2;
198 
199  if (nSigmaCut == 0) break; // no rejection step to perform
200  MeasuredStarList msOutliers;
201  FittedStarList fsOutliers;
202  // keep nOutliers so we don't have to sum msOutliers.size()+fsOutliers.size() twice below.
203  int nOutliers = findOutliers(nSigmaCut, msOutliers, fsOutliers);
204  totalMeasOutliers += msOutliers.size();
205  totalRefOutliers += fsOutliers.size();
206  if (nOutliers == 0) break;
207  TripletList outlierTriplets(nOutliers);
208  grad.setZero(); // recycle the gradient
209  // compute the contributions of outliers to derivatives
210  outliersContributions(msOutliers, fsOutliers, outlierTriplets, grad);
211  // Remove significant outliers
212  removeMeasOutliers(msOutliers);
213  removeRefOutliers(fsOutliers);
214  if (doRankUpdate) {
215  // convert triplet list to eigen internal format
216  SparseMatrixD H(_nParTot, outlierTriplets.getNextFreeIndex());
217  H.setFromTriplets(outlierTriplets.begin(), outlierTriplets.end());
218  chol.update(H, false /* means downdate */);
219  // The contribution of outliers to the gradient is the opposite
220  // of the contribution of all other terms, because they add up to 0
221  grad *= -1;
222  } else {
223  // don't reuse tripletList because we want a new nextFreeIndex.
224  TripletList nextTripletList(_lastNTrip);
225  grad.setZero();
226  // Rebuild the matrix and gradient
227  leastSquareDerivatives(nextTripletList, grad);
228  _lastNTrip = nextTripletList.size();
229  LOGLS_DEBUG(_log, "Triplets recomputed, ntrip = " << nextTripletList.size());
230 
231  hessian = createHessian(_nParTot, nextTripletList);
232  nextTripletList.clear(); // we don't need it any more after we have the hessian.
233 
235  "Restarting factorization, hessian: dim="
236  << hessian.rows() << " non-zeros=" << hessian.nonZeros()
237  << " filling-frac = " << hessian.nonZeros() / std::pow(hessian.rows(), 2));
238  chol.compute(hessian);
239  if (chol.info() != Eigen::Success) {
240  LOGLS_ERROR(_log, "minimize: factorization failed ");
241  return MinimizeResult::Failed;
242  }
243  }
244  }
245 
246  // only print the outlier summary if outlier rejection was turned on.
247  if (nSigmaCut != 0) {
248  LOGLS_INFO(_log, "Number of outliers (Measured + Reference = Total): "
249  << totalMeasOutliers << " + " << totalRefOutliers << " = "
250  << totalMeasOutliers + totalRefOutliers);
251  }
252  return returnCode;
253 }
254 
256  TripletList &tripletList, Eigen::VectorXd &grad) {
257  for (auto &outlier : msOutliers) {
258  MeasuredStarList tmp;
259  tmp.push_back(outlier);
260  const CcdImage &ccdImage = outlier->getCcdImage();
261  leastSquareDerivativesMeasurement(ccdImage, tripletList, grad, &tmp);
262  }
263  leastSquareDerivativesReference(fsOutliers, tripletList, grad);
264 }
265 
267  for (auto &measuredStar : outliers) {
268  auto fittedStar = measuredStar->getFittedStar();
269  measuredStar->setValid(false);
270  fittedStar->getMeasurementCount()--; // could be put in setValid
271  }
272 }
273 
275  for (auto &fittedStar : outliers) {
276  fittedStar->setRefStar(nullptr);
277  }
278 }
279 
280 void FitterBase::leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const {
281  auto ccdImageList = _associations->getCcdImageList();
282  for (auto const &ccdImage : ccdImageList) {
283  leastSquareDerivativesMeasurement(*ccdImage, tripletList, grad);
284  }
285  leastSquareDerivativesReference(_associations->fittedStarList, tripletList, grad);
286 }
287 
288 void FitterBase::saveChi2Contributions(std::string const &baseName) const {
289  /* cook-up 2 different file names by inserting something just before
290  the dot (if any), and within the actual file name. */
291  size_t dot = baseName.rfind('.');
292  size_t slash = baseName.rfind('/');
293  if (dot == std::string::npos || (slash != std::string::npos && dot < slash)) dot = baseName.size();
294  std::string measTuple(baseName);
295  measTuple.insert(dot, "-meas");
296  saveChi2MeasContributions(measTuple);
297  std::string refTuple(baseName);
298  refTuple.insert(dot, "-ref");
299  saveChi2RefContributions(refTuple);
300 }
301 
302 } // namespace jointcal
303 } // namespace lsst
virtual void accumulateStatRefStars(Chi2Accumulator &accum) const =0
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars...
virtual void leastSquareDerivativesReference(FittedStarList const &fittedStarList, TripletList &tripletList, Eigen::VectorXd &grad) const =0
Compute the derivatives of the reference terms.
Simple structure to accumulate chi2 and ndof.
Definition: Chi2.h:29
Eigen::SparseMatrix< double > SparseMatrixD
Definition: Eigenstuff.h:12
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...
unsigned findOutliers(double nSigmaCut, MeasuredStarList &msOutliers, FittedStarList &fsOutliers) const
Find Measurements and references contributing more than a cut, computed as The outliers are NOT remo...
Definition: FitterBase.cc:26
T rend(T... args)
T rfind(T... args)
T endl(T... args)
virtual void getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, std::vector< unsigned > &indices) const =0
Set the indices of a measured star from the full matrix, for outlier removal.
T end(T... args)
void removeMeasOutliers(MeasuredStarList &outliers)
Remove measuredStar outliers from the fit. No Refit done.
Definition: FitterBase.cc:266
void update(SparseMatrixD const &H, bool UpOrDown)
Definition: Eigenstuff.h:42
A list of MeasuredStar. They are usually filled in Associations::createCcdImage.
Definition: MeasuredStar.h:113
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:288
MinimizeResult
Return value of minimize()
Definition: FitterBase.h:17
#define LOGLS_INFO(logger, message)
def dot(symb, c, r, frame=None, size=2, ctype=None, origin=afwImage.PARENT, args, kwargs)
STL class.
virtual void saveChi2MeasContributions(std::string const &baseName) const =0
Save a CSV file containing residuals of measurement terms.
Structure to accumulate the chi2 contributions per each star (to help find outliers).
Definition: Chi2.h:77
T push_back(T... args)
STL class.
Class for a simple mapping implementing a generic Gtransfo.
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition: FittedStar.h:121
MinimizeResult minimize(std::string const &whatToFit, double const nSigmaCut=0, bool const doRankUpdate=true, std::string const &dumpMatrixFile="")
Does a 1 step minimization, assuming a linear model.
Definition: FitterBase.cc:143
void removeRefOutliers(FittedStarList &outliers)
Remove refStar outliers from the fit. No Refit done.
Definition: FitterBase.cc:274
objects measured on actual images.
Definition: MeasuredStar.h:19
T dynamic_pointer_cast(T... args)
T clear(T... args)
virtual void saveChi2RefContributions(std::string const &baseName) const =0
Save a CSV file containing residuals of reference terms.
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:129
virtual void assignIndices(std::string const &whatToFit)=0
Set parameters to fit and assign indices in the big matrix.
T insert(T... args)
T size(T... args)
Chi2Statistic computeChi2() const
Returns the chi2 for the current state.
Definition: FitterBase.cc:16
T begin(T... args)
T pow(T... args)
#define LOGLS_DEBUG(logger, message)
void leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const
Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting...
Definition: FitterBase.cc:280
T sort(T... args)
Handler of an actual image from a single CCD.
Definition: CcdImage.h:35
#define LOGL_WARN(logger, message...)
virtual void offsetParams(Eigen::VectorXd const &delta)=0
Offset the parameters by the requested quantities.
unsigned getNextFreeIndex() const
Definition: Tripletlist.h:22
void outliersContributions(MeasuredStarList &msOutliers, FittedStarList &fsOutliers, TripletList &tripletList, Eigen::VectorXd &grad)
Contributions to derivatives from (presumably) outlier terms.
Definition: FitterBase.cc:255
The objects which have been measured several times.
Definition: FittedStar.h:37
std::pair< double, double > computeAverageAndSigma()
Compute the average and std-deviation of these chisq values.
Definition: Chi2.cc:31
#define LOGLS_ERROR(logger, message)
#define LOGLS_WARN(logger, message)
T reserve(T... args)
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.
T rbegin(T... args)