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