lsst.jointcal  14.0-20-g39ceace+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 MinimizeResult FitterBase::minimize(std::string const &whatToFit, double nSigmaCut) {
120  assignIndices(whatToFit);
121 
123 
124  // TODO : write a guesser for the number of triplets
125  unsigned nTrip = (_lastNTrip) ? _lastNTrip : 1e6;
126  TripletList tripletList(nTrip);
127  Eigen::VectorXd grad(_nParTot);
128  grad.setZero();
129 
130  // Fill the triplets
131  leastSquareDerivatives(tripletList, grad);
132  _lastNTrip = tripletList.size();
133 
134  LOGLS_DEBUG(_log, "End of triplet filling, ntrip = " << tripletList.size());
135 
136  SpMat hessian;
137  {
138  SpMat jacobian(_nParTot, tripletList.getNextFreeIndex());
139  jacobian.setFromTriplets(tripletList.begin(), tripletList.end());
140  // release memory shrink_to_fit is C++11
141  tripletList.clear(); // tripletList.shrink_to_fit();
142  hessian = jacobian * jacobian.transpose();
143  } // release the Jacobian
144 
145  LOGLS_DEBUG(_log, "Starting factorization, hessian: dim="
146  << hessian.rows() << " non-zeros=" << hessian.nonZeros()
147  << " filling-frac = " << hessian.nonZeros() / std::pow(hessian.rows(), 2));
148 
149  CholmodSimplicialLDLT2<SpMat> chol(hessian);
150  if (chol.info() != Eigen::Success) {
151  LOGLS_ERROR(_log, "minimize: factorization failed ");
152  return MinimizeResult::Failed;
153  }
154 
155  unsigned totalOutliers = 0;
156  double oldChi2 = computeChi2().chi2;
157 
158  while (true) {
159  Eigen::VectorXd delta = chol.solve(grad);
160  offsetParams(delta);
161  Chi2Statistic currentChi2(computeChi2());
162  LOGLS_DEBUG(_log, currentChi2);
163  if (currentChi2.chi2 > oldChi2 && totalOutliers != 0) {
164  LOGL_WARN(_log, "chi2 went up, skipping outlier rejection loop");
165  returnCode = MinimizeResult::Chi2Increased;
166  break;
167  }
168  oldChi2 = currentChi2.chi2;
169 
170  if (nSigmaCut == 0) break; // no rejection step to perform
171  MeasuredStarList msOutliers;
172  FittedStarList fsOutliers;
173  int nOutliers = findOutliers(nSigmaCut, msOutliers, fsOutliers);
174  totalOutliers += nOutliers;
175  if (nOutliers == 0) break;
176  TripletList tripletList(nOutliers);
177  grad.setZero(); // recycle the gradient
178  // compute the contributions of outliers to derivatives
179  outliersContributions(msOutliers, fsOutliers, tripletList, grad);
180  // Remove significant outliers
181  removeMeasOutliers(msOutliers);
182  removeRefOutliers(fsOutliers);
183  // convert triplet list to eigen internal format
184  SpMat H(_nParTot, tripletList.getNextFreeIndex());
185  H.setFromTriplets(tripletList.begin(), tripletList.end());
186  int update_status = chol.update(H, false /* means downdate */);
187  LOGLS_DEBUG(_log, "cholmod update_status " << update_status);
188  // The contribution of outliers to the gradient is the opposite
189  // of the contribution of all other terms, because they add up to 0
190  grad *= -1;
191  }
192 
193  LOGLS_INFO(_log, "Total number of outliers " << totalOutliers);
194  return returnCode;
195 }
196 
198  TripletList &tripletList, Eigen::VectorXd &grad) {
199  for (auto &outlier : msOutliers) {
200  MeasuredStarList tmp;
201  tmp.push_back(outlier);
202  const CcdImage &ccdImage = outlier->getCcdImage();
203  leastSquareDerivativesMeasurement(ccdImage, tripletList, grad, &tmp);
204  }
205  leastSquareDerivativesReference(fsOutliers, tripletList, grad);
206 }
207 
209  for (auto &measuredStar : outliers) {
210  auto fittedStar = std::const_pointer_cast<FittedStar>(measuredStar->getFittedStar());
211  measuredStar->setValid(false);
212  fittedStar->getMeasurementCount()--; // could be put in setValid
213  }
214 }
215 
217  for (auto &fittedStar : outliers) {
218  fittedStar->setRefStar(nullptr);
219  }
220 }
221 
222 void FitterBase::leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const {
223  auto ccdImageList = _associations->getCcdImageList();
224  for (auto const &ccdImage : ccdImageList) {
225  leastSquareDerivativesMeasurement(*ccdImage, tripletList, grad);
226  }
227  leastSquareDerivativesReference(_associations->fittedStarList, tripletList, grad);
228 }
229 
230 void FitterBase::saveChi2Contributions(std::string const &baseName) const {
231  /* cook-up 2 different file names by inserting something just before
232  the dot (if any), and within the actual file name. */
233  size_t dot = baseName.rfind('.');
234  size_t slash = baseName.rfind('/');
235  if (dot == std::string::npos || (slash != std::string::npos && dot < slash)) dot = baseName.size();
236  std::string measTuple(baseName);
237  measTuple.insert(dot, "-meas");
238  saveChi2MeasContributions(measTuple);
239  std::string refTuple(baseName);
240  refTuple.insert(dot, "-ref");
241  saveChi2RefContributions(refTuple);
242 }
243 
244 } // namespace jointcal
245 } // namespace lsst
#define LOGLS_WARN(logger, message)
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
int update(SpMat const &H, bool UpOrDown)
Definition: Eigenstuff.h:41
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)
int getMeasurementCount() const
Definition: FittedStar.h:82
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.
MinimizeResult minimize(std::string const &whatToFit, double nSigmaCut=0)
Does a 1 step minimization, assuming a linear model.
Definition: FitterBase.cc:119
T end(T... args)
void removeMeasOutliers(MeasuredStarList &outliers)
Remove measuredStar outliers from the fit. No Refit done.
Definition: FitterBase.cc:208
A list of MeasuredStar. They are usually filled in Associations::AddImage.
Definition: MeasuredStar.h:111
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:230
MinimizeResult
Return value of minimize()
Definition: FitterBase.h:17
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)
#define LOGLS_DEBUG(logger, message)
Class for a simple mapping implementing a generic Gtransfo.
Eigen::SparseMatrix< double > SpMat
Definition: Eigenstuff.h:11
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:216
objects measured on actual images.
Definition: MeasuredStar.h:18
T dynamic_pointer_cast(T... args)
#define LOGL_WARN(logger, message...)
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:114
virtual void assignIndices(std::string const &whatToFit)=0
Set parameters to fit and assign indices in the big matrix.
#define LOGLS_INFO(logger, message)
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)
void leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const
Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting...
Definition: FitterBase.cc:222
T sort(T... args)
Handler of an actual image from a single CCD.
Definition: CcdImage.h:35
virtual void offsetParams(Eigen::VectorXd const &delta)=0
Offset the parameters by the requested quantities.
#define LOGLS_ERROR(logger, message)
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:197
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
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)