LSSTApplications  18.0.0+106,18.0.0+50,19.0.0,19.0.0+1,19.0.0+10,19.0.0+11,19.0.0+13,19.0.0+17,19.0.0+2,19.0.0-1-g20d9b18+6,19.0.0-1-g425ff20,19.0.0-1-g5549ca4,19.0.0-1-g580fafe+6,19.0.0-1-g6fe20d0+1,19.0.0-1-g7011481+9,19.0.0-1-g8c57eb9+6,19.0.0-1-gb5175dc+11,19.0.0-1-gdc0e4a7+9,19.0.0-1-ge272bc4+6,19.0.0-1-ge3aa853,19.0.0-10-g448f008b,19.0.0-12-g6990b2c,19.0.0-2-g0d9f9cd+11,19.0.0-2-g3d9e4fb2+11,19.0.0-2-g5037de4,19.0.0-2-gb96a1c4+3,19.0.0-2-gd955cfd+15,19.0.0-3-g2d13df8,19.0.0-3-g6f3c7dc,19.0.0-4-g725f80e+11,19.0.0-4-ga671dab3b+1,19.0.0-4-gad373c5+3,19.0.0-5-ga2acb9c+2,19.0.0-5-gfe96e6c+2,w.2020.01
LSSTDataManagementBasePackage
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
lsst::jointcal::FitterBase Class Referenceabstract

Base class for fitters. More...

#include <FitterBase.h>

Inheritance diagram for lsst::jointcal::FitterBase:
lsst::jointcal::AstrometryFit lsst::jointcal::PhotometryFit

Public Member Functions

 FitterBase (std::shared_ptr< Associations > associations)
 
 FitterBase (FitterBase const &)=delete
 No copy or move: there is only ever one fitter of a given type. More...
 
 FitterBase (FitterBase &&)=delete
 
FitterBaseoperator= (FitterBase const &)=delete
 
FitterBaseoperator= (FitterBase &&)=delete
 
MinimizeResult minimize (std::string const &whatToFit, double const nSigmaCut=0, bool const doRankUpdate=true, bool const doLineSearch=false, std::string const &dumpMatrixFile="")
 Does a 1 step minimization, assuming a linear model. More...
 
Chi2Statistic computeChi2 () const
 Returns the chi2 for the current state. More...
 
void leastSquareDerivatives (TripletList &tripletList, Eigen::VectorXd &grad) const
 Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting. More...
 
virtual void offsetParams (Eigen::VectorXd const &delta)=0
 Offset the parameters by the requested quantities. More...
 
virtual void assignIndices (std::string const &whatToFit)=0
 Set parameters to fit and assign indices in the big matrix. More...
 
virtual void saveChi2Contributions (std::string const &baseName) const
 Save the full chi2 term per star that was used in the minimization, for debugging. More...
 

Protected Member Functions

virtual void saveChi2MeasContributions (std::string const &filename) const =0
 Save a CSV file containing residuals of measurement terms. More...
 
virtual void saveChi2RefContributions (std::string const &filename) const =0
 Save a CSV file containing residuals of reference terms. More...
 
std::size_t findOutliers (double nSigmaCut, MeasuredStarList &msOutliers, FittedStarList &fsOutliers) const
 Find Measurements and references contributing more than a cut, computed as

\[ <chi2> + nSigmaCut + rms(chi2). \]

The outliers are NOT removed, and no refit is done. More...

 
void outliersContributions (MeasuredStarList &msOutliers, FittedStarList &fsOutliers, TripletList &tripletList, Eigen::VectorXd &grad)
 Contributions to derivatives from (presumably) outlier terms. More...
 
void removeMeasOutliers (MeasuredStarList &outliers)
 Remove measuredStar outliers from the fit. No Refit done. More...
 
void removeRefOutliers (FittedStarList &outliers)
 Remove refStar outliers from the fit. No Refit done. More...
 
virtual void getIndicesOfMeasuredStar (MeasuredStar const &measuredStar, IndexVector &indices) const =0
 Set the indices of a measured star from the full matrix, for outlier removal. More...
 
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. More...
 
virtual void accumulateStatRefStars (Chi2Accumulator &accum) const =0
 Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars. More...
 
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. More...
 
virtual void leastSquareDerivativesReference (FittedStarList const &fittedStarList, TripletList &tripletList, Eigen::VectorXd &grad) const =0
 Compute the derivatives of the reference terms. More...
 

Protected Attributes

std::shared_ptr< Associations_associations
 
std::string _whatToFit
 
Eigen::Index _lastNTrip
 
Eigen::Index _nParTot
 
Eigen::Index _nMeasuredStars
 
LOG_LOGGER _log
 

Detailed Description

Base class for fitters.

Implements minimize and findOutliers. Chi2, residual, derivative, etc. calculations must be implemented in the child class via the virtual methods.

Definition at line 53 of file FitterBase.h.

Constructor & Destructor Documentation

◆ FitterBase() [1/3]

lsst::jointcal::FitterBase::FitterBase ( std::shared_ptr< Associations associations)
inlineexplicit

Definition at line 55 of file FitterBase.h.

56  : _associations(associations), _whatToFit(""), _lastNTrip(0), _nParTot(0), _nMeasuredStars(0) {}
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:153
Eigen::Index _nMeasuredStars
Definition: FitterBase.h:158

◆ FitterBase() [2/3]

lsst::jointcal::FitterBase::FitterBase ( FitterBase const &  )
delete

No copy or move: there is only ever one fitter of a given type.

◆ FitterBase() [3/3]

lsst::jointcal::FitterBase::FitterBase ( FitterBase &&  )
delete

Member Function Documentation

◆ accumulateStatImageList()

virtual void lsst::jointcal::FitterBase::accumulateStatImageList ( CcdImageList const &  ccdImageList,
Chi2Accumulator accum 
) const
protectedpure virtual

Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for measurements.

◆ accumulateStatRefStars()

virtual void lsst::jointcal::FitterBase::accumulateStatRefStars ( Chi2Accumulator accum) const
protectedpure virtual

Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars.

◆ assignIndices()

virtual void lsst::jointcal::FitterBase::assignIndices ( std::string const &  whatToFit)
pure virtual

Set parameters to fit and assign indices in the big matrix.

Parameters
[in]whatToFitSee child class documentation for valid string values.

Implemented in lsst::jointcal::AstrometryFit, and lsst::jointcal::PhotometryFit.

◆ computeChi2()

Chi2Statistic lsst::jointcal::FitterBase::computeChi2 ( ) const

Returns the chi2 for the current state.

Definition at line 42 of file FitterBase.cc.

42  {
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 }
virtual void accumulateStatRefStars(Chi2Accumulator &accum) const =0
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars...
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...
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:153

◆ findOutliers()

std::size_t lsst::jointcal::FitterBase::findOutliers ( double  nSigmaCut,
MeasuredStarList msOutliers,
FittedStarList fsOutliers 
) const
protected

Find Measurements and references contributing more than a cut, computed as

\[ <chi2> + nSigmaCut + rms(chi2). \]

The outliers are NOT removed, and no refit is done.

After returning from here, there are still measurements that contribute above the cut, but their contribution should be evaluated after a refit before discarding them.

Parameters
[in]nSigmaCutNumber of sigma to select on.
[out]msOutlierslist of MeasuredStar outliers to populate
[out]fsOutlierslist of FittedStar outliers to populate
Returns
Total number of outliers that were removed.

Definition at line 52 of file FitterBase.cc.

53  {
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  double 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 }
#define LOGLS_WARN(logger, message)
Log a warn-level message using an iostream-based interface.
Definition: Log.h:648
virtual void accumulateStatRefStars(Chi2Accumulator &accum) const =0
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars...
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...
T push_back(T... args)
#define LOGLS_DEBUG(logger, message)
Log a debug-level message using an iostream-based interface.
Definition: Log.h:608
#define LOGLS_TRACE(logger, message)
Log a trace-level message using an iostream-based interface.
Definition: Log.h:588
T dynamic_pointer_cast(T... args)
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:153
#define LOGLS_INFO(logger, message)
Log a info-level message using an iostream-based interface.
Definition: Log.h:628
STL class.
T sort(T... args)
Eigen::Index _nMeasuredStars
Definition: FitterBase.h:158
virtual void getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, IndexVector &indices) const =0
Set the indices of a measured star from the full matrix, for outlier removal.

◆ getIndicesOfMeasuredStar()

virtual void lsst::jointcal::FitterBase::getIndicesOfMeasuredStar ( MeasuredStar const &  measuredStar,
IndexVector indices 
) const
protectedpure virtual

Set the indices of a measured star from the full matrix, for outlier removal.

◆ leastSquareDerivatives()

void lsst::jointcal::FitterBase::leastSquareDerivatives ( TripletList tripletList,
Eigen::VectorXd &  grad 
) const

Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting.

The Jacobian is given as triplets in a sparse matrix, the gradient as a dense vector. The parameters which vary, and their indices, are to be set using assignIndices.

Parameters
tripletListtripletList of (row,col,value) representing the Jacobian of the chi2.
gradThe gradient of the chi2.

Definition at line 317 of file FitterBase.cc.

317  {
318  auto ccdImageList = _associations->getCcdImageList();
319  for (auto const &ccdImage : ccdImageList) {
320  leastSquareDerivativesMeasurement(*ccdImage, tripletList, grad);
321  }
322  leastSquareDerivativesReference(_associations->fittedStarList, tripletList, grad);
323 }
virtual void leastSquareDerivativesReference(FittedStarList const &fittedStarList, TripletList &tripletList, Eigen::VectorXd &grad) const =0
Compute the derivatives of the reference terms.
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:153
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.

◆ leastSquareDerivativesMeasurement()

virtual void lsst::jointcal::FitterBase::leastSquareDerivativesMeasurement ( CcdImage const &  ccdImage,
TripletList tripletList,
Eigen::VectorXd &  grad,
MeasuredStarList const *  measuredStarList = nullptr 
) const
protectedpure virtual

Compute the derivatives of the measured stars and model for one CcdImage.

The last argument will process a sub-list for outlier removal.

◆ leastSquareDerivativesReference()

virtual void lsst::jointcal::FitterBase::leastSquareDerivativesReference ( FittedStarList const &  fittedStarList,
TripletList tripletList,
Eigen::VectorXd &  grad 
) const
protectedpure virtual

Compute the derivatives of the reference terms.

◆ minimize()

MinimizeResult lsst::jointcal::FitterBase::minimize ( std::string const &  whatToFit,
double const  nSigmaCut = 0,
bool const  doRankUpdate = true,
bool const  doLineSearch = false,
std::string const &  dumpMatrixFile = "" 
)

Does a 1 step minimization, assuming a linear model.

This is a complete Newton Raphson step. Compute first and second derivatives, solve for the step and apply it, with an optional line search.

It calls assignIndices, leastSquareDerivatives, solves the linear system and calls offsetParams, then removes outliers in a loop if requested. Relies on sparse linear algebra via Eigen's CholmodSupport package.

Parameters
[in]whatToFitSee child method assignIndices for valid string values.
[in]nSigmaCutHow many sigma to reject outliers at. Outlier rejection ignored for nSigmaCut=0.
[in]doRankUpdateUse CholmodSimplicialLDLT2.update() to do a fast rank update after outlier removal; otherwise do a slower full recomputation of the matrix. Only matters if nSigmaCut != 0.
[in]doLineSearchUse boost's brent_find_minima to perform a line search after the gradient solution is found, and apply the scale factor to the computed offsets. The line search is done in the domain [-1, 2], but if the scale factor is far from 1.0, then the problem is likely in a significantly non-linear regime.
[in]dumpMatrixFileWrite the pre-fit Hessian matrix and gradient to the files with "-mat.txt" and "-grad.txt". Be aware, this requires a large increase in memory usage to create a dense matrix before writing it; the output file may be large. Writing the matrix can be helpful for debugging bad fits. Read it and compute the real eigenvalues (recall that the Hessian is symmetric by construction) with numpy:
hessian = np.matrix(np.loadtxt("dumpMatrixFile-mat.txt"))
values, vectors = np.linalg.eigh(hessian)
Returns
Return code describing success/failure of fit.
Note
When fitting one parameter set by itself (e.g. "Model"), the system is purely linear (assuming there are no cross-terms in the derivatives, e.g. the SimpleAstrometryModel), which should result in the optimal chi2 after a single step. This can be used to debug the fitter by fitting that parameter set twice in a row: the second run with the same "whatToFit" will produce no change in the fitted parameters, if the calculations and indices are defined correctly.

Definition at line 171 of file FitterBase.cc.

172  {
173  assignIndices(whatToFit);
174 
176 
177  // TODO : write a guesser for the number of triplets
178  std::size_t nTrip = (_lastNTrip) ? _lastNTrip : 1e6;
179  TripletList tripletList(nTrip);
180  Eigen::VectorXd grad(_nParTot);
181  grad.setZero();
182  double scale = 1.0;
183 
184  // Fill the triplets
185  leastSquareDerivatives(tripletList, grad);
186  _lastNTrip = tripletList.size();
187 
188  LOGLS_DEBUG(_log, "End of triplet filling, ntrip = " << tripletList.size());
189 
190  SparseMatrixD hessian = createHessian(_nParTot, tripletList);
191  tripletList.clear(); // we don't need it any more after we have the hessian.
192 
193  LOGLS_DEBUG(_log, "Starting factorization, hessian: dim="
194  << hessian.rows() << " non-zeros=" << hessian.nonZeros()
195  << " filling-frac = " << hessian.nonZeros() / std::pow(hessian.rows(), 2));
196 
197  if (dumpMatrixFile != "") {
198  if (hessian.rows() * hessian.cols() > 2e8) {
199  LOGLS_WARN(_log, "Hessian matrix is too big to dump to file, with rows, columns: "
200  << hessian.rows() << ", " << hessian.cols());
201  } else {
202  dumpMatrixAndGradient(hessian, grad, dumpMatrixFile, _log);
203  }
204  }
205 
207  if (chol.info() != Eigen::Success) {
208  LOGLS_ERROR(_log, "minimize: factorization failed ");
209  return MinimizeResult::Failed;
210  }
211 
212  std::size_t totalMeasOutliers = 0;
213  std::size_t totalRefOutliers = 0;
214  double oldChi2 = computeChi2().chi2;
215 
216  while (true) {
217  Eigen::VectorXd delta = chol.solve(grad);
218  if (doLineSearch) {
219  scale = _lineSearch(delta);
220  }
221  offsetParams(scale * delta);
222  Chi2Statistic currentChi2(computeChi2());
223  LOGLS_DEBUG(_log, currentChi2);
224  if (!isfinite(currentChi2.chi2)) {
225  LOGL_ERROR(_log, "chi2 is not finite. Aborting outlier rejection.");
226  returnCode = MinimizeResult::NonFinite;
227  break;
228  }
229  if (currentChi2.chi2 > oldChi2 && totalMeasOutliers + totalRefOutliers != 0) {
230  LOGL_WARN(_log, "chi2 went up, skipping outlier rejection loop");
231  returnCode = MinimizeResult::Chi2Increased;
232  break;
233  }
234  oldChi2 = currentChi2.chi2;
235 
236  if (nSigmaCut == 0) break; // no rejection step to perform
237  MeasuredStarList msOutliers;
238  FittedStarList fsOutliers;
239  // keep nOutliers so we don't have to sum msOutliers.size()+fsOutliers.size() twice below.
240  std::size_t nOutliers = findOutliers(nSigmaCut, msOutliers, fsOutliers);
241  totalMeasOutliers += msOutliers.size();
242  totalRefOutliers += fsOutliers.size();
243  if (nOutliers == 0) break;
244  TripletList outlierTriplets(nOutliers);
245  grad.setZero(); // recycle the gradient
246  // compute the contributions of outliers to derivatives
247  outliersContributions(msOutliers, fsOutliers, outlierTriplets, grad);
248  // Remove significant outliers
249  removeMeasOutliers(msOutliers);
250  removeRefOutliers(fsOutliers);
251  if (doRankUpdate) {
252  // convert triplet list to eigen internal format
253  SparseMatrixD H(_nParTot, outlierTriplets.getNextFreeIndex());
254  H.setFromTriplets(outlierTriplets.begin(), outlierTriplets.end());
255  chol.update(H, false /* means downdate */);
256  // The contribution of outliers to the gradient is the opposite
257  // of the contribution of all other terms, because they add up to 0
258  grad *= -1;
259  } else {
260  // don't reuse tripletList because we want a new nextFreeIndex.
261  TripletList nextTripletList(_lastNTrip);
262  grad.setZero();
263  // Rebuild the matrix and gradient
264  leastSquareDerivatives(nextTripletList, grad);
265  _lastNTrip = nextTripletList.size();
266  LOGLS_DEBUG(_log, "Triplets recomputed, ntrip = " << nextTripletList.size());
267 
268  hessian = createHessian(_nParTot, nextTripletList);
269  nextTripletList.clear(); // we don't need it any more after we have the hessian.
270 
272  "Restarting factorization, hessian: dim="
273  << hessian.rows() << " non-zeros=" << hessian.nonZeros()
274  << " filling-frac = " << hessian.nonZeros() / std::pow(hessian.rows(), 2));
275  chol.compute(hessian);
276  if (chol.info() != Eigen::Success) {
277  LOGLS_ERROR(_log, "minimize: factorization failed ");
278  return MinimizeResult::Failed;
279  }
280  }
281  }
282 
283  // only print the outlier summary if outlier rejection was turned on.
284  if (nSigmaCut != 0) {
285  LOGLS_INFO(_log, "Number of outliers (Measured + Reference = Total): "
286  << totalMeasOutliers << " + " << totalRefOutliers << " = "
287  << totalMeasOutliers + totalRefOutliers);
288  }
289  return returnCode;
290 }
#define LOGLS_WARN(logger, message)
Log a warn-level message using an iostream-based interface.
Definition: Log.h:648
std::size_t 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:52
def scale(algorithm, min, max=None, frame=None)
Definition: ds9.py:109
#define LOGL_ERROR(logger, message...)
Log a error-level message using a varargs/printf style interface.
Definition: Log.h:552
void removeMeasOutliers(MeasuredStarList &outliers)
Remove measuredStar outliers from the fit. No Refit done.
Definition: FitterBase.cc:303
MinimizeResult
Return value of minimize()
Definition: FitterBase.h:40
Eigen::SparseMatrix< double, 0, Eigen::Index > SparseMatrixD
Definition: Eigenstuff.h:35
#define LOGLS_DEBUG(logger, message)
Log a debug-level message using an iostream-based interface.
Definition: Log.h:608
void removeRefOutliers(FittedStarList &outliers)
Remove refStar outliers from the fit. No Refit done.
Definition: FitterBase.cc:311
T isfinite(T... args)
#define LOGL_WARN(logger, message...)
Log a warn-level message using a varargs/printf style interface.
Definition: Log.h:536
virtual void assignIndices(std::string const &whatToFit)=0
Set parameters to fit and assign indices in the big matrix.
#define LOGLS_INFO(logger, message)
Log a info-level message using an iostream-based interface.
Definition: Log.h:628
Chi2Statistic computeChi2() const
Returns the chi2 for the current state.
Definition: FitterBase.cc:42
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:317
virtual void offsetParams(Eigen::VectorXd const &delta)=0
Offset the parameters by the requested quantities.
#define LOGLS_ERROR(logger, message)
Log a error-level message using an iostream-based interface.
Definition: Log.h:668
void outliersContributions(MeasuredStarList &msOutliers, FittedStarList &fsOutliers, TripletList &tripletList, Eigen::VectorXd &grad)
Contributions to derivatives from (presumably) outlier terms.
Definition: FitterBase.cc:292

◆ offsetParams()

virtual void lsst::jointcal::FitterBase::offsetParams ( Eigen::VectorXd const &  delta)
pure virtual

Offset the parameters by the requested quantities.

The used parameter layout is the one from the last call to assignIndices or minimize(). There is no easy way to check that the current setting of whatToFit and the provided Delta vector are compatible: we can only test the size.

Parameters
[in]deltavector of offsets to apply

Implemented in lsst::jointcal::AstrometryFit, and lsst::jointcal::PhotometryFit.

◆ operator=() [1/2]

FitterBase& lsst::jointcal::FitterBase::operator= ( FitterBase const &  )
delete

◆ operator=() [2/2]

FitterBase& lsst::jointcal::FitterBase::operator= ( FitterBase &&  )
delete

◆ outliersContributions()

void lsst::jointcal::FitterBase::outliersContributions ( MeasuredStarList msOutliers,
FittedStarList fsOutliers,
TripletList tripletList,
Eigen::VectorXd &  grad 
)
protected

Contributions to derivatives from (presumably) outlier terms.

No discarding done.

Definition at line 292 of file FitterBase.cc.

293  {
294  for (auto &outlier : msOutliers) {
295  MeasuredStarList tmp;
296  tmp.push_back(outlier);
297  const CcdImage &ccdImage = outlier->getCcdImage();
298  leastSquareDerivativesMeasurement(ccdImage, tripletList, grad, &tmp);
299  }
300  leastSquareDerivativesReference(fsOutliers, tripletList, grad);
301 }
virtual void leastSquareDerivativesReference(FittedStarList const &fittedStarList, TripletList &tripletList, Eigen::VectorXd &grad) const =0
Compute the derivatives of the 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.

◆ removeMeasOutliers()

void lsst::jointcal::FitterBase::removeMeasOutliers ( MeasuredStarList outliers)
protected

Remove measuredStar outliers from the fit. No Refit done.

Definition at line 303 of file FitterBase.cc.

303  {
304  for (auto &measuredStar : outliers) {
305  auto fittedStar = measuredStar->getFittedStar();
306  measuredStar->setValid(false);
307  fittedStar->getMeasurementCount()--; // could be put in setValid
308  }
309 }

◆ removeRefOutliers()

void lsst::jointcal::FitterBase::removeRefOutliers ( FittedStarList outliers)
protected

Remove refStar outliers from the fit. No Refit done.

Definition at line 311 of file FitterBase.cc.

311  {
312  for (auto &fittedStar : outliers) {
313  fittedStar->setRefStar(nullptr);
314  }
315 }

◆ saveChi2Contributions()

void lsst::jointcal::FitterBase::saveChi2Contributions ( std::string const &  baseName) const
virtual

Save the full chi2 term per star that was used in the minimization, for debugging.

Saves results to text files "baseName-meas.csv" and "baseName-ref.csv" for the MeasuredStar and RefStar contributions, respectively. This method is mostly useful for debugging: we will probably want to create a better persistence system for jointcal's internal representations in the future (see DM-12446).

Definition at line 325 of file FitterBase.cc.

325  {
326  std::string replaceStr = "{type}";
327  auto pos = baseName.find(replaceStr);
328  std::string measFilename(baseName);
329  measFilename.replace(pos, replaceStr.size(), "-meas.csv");
330  std::string refFilename(baseName);
331  refFilename.replace(pos, replaceStr.size(), "-ref.csv");
332  saveChi2MeasContributions(measFilename);
333  saveChi2RefContributions(refFilename);
334 }
virtual void saveChi2RefContributions(std::string const &filename) const =0
Save a CSV file containing residuals of reference terms.
STL class.
virtual void saveChi2MeasContributions(std::string const &filename) const =0
Save a CSV file containing residuals of measurement terms.
T find(T... args)
T size(T... args)

◆ saveChi2MeasContributions()

virtual void lsst::jointcal::FitterBase::saveChi2MeasContributions ( std::string const &  filename) const
protectedpure virtual

Save a CSV file containing residuals of measurement terms.

Implemented in lsst::jointcal::AstrometryFit, and lsst::jointcal::PhotometryFit.

◆ saveChi2RefContributions()

virtual void lsst::jointcal::FitterBase::saveChi2RefContributions ( std::string const &  filename) const
protectedpure virtual

Save a CSV file containing residuals of reference terms.

Implemented in lsst::jointcal::AstrometryFit, and lsst::jointcal::PhotometryFit.

Member Data Documentation

◆ _associations

std::shared_ptr<Associations> lsst::jointcal::FitterBase::_associations
protected

Definition at line 153 of file FitterBase.h.

◆ _lastNTrip

Eigen::Index lsst::jointcal::FitterBase::_lastNTrip
protected

Definition at line 156 of file FitterBase.h.

◆ _log

LOG_LOGGER lsst::jointcal::FitterBase::_log
protected

Definition at line 161 of file FitterBase.h.

◆ _nMeasuredStars

Eigen::Index lsst::jointcal::FitterBase::_nMeasuredStars
protected

Definition at line 158 of file FitterBase.h.

◆ _nParTot

Eigen::Index lsst::jointcal::FitterBase::_nParTot
protected

Definition at line 157 of file FitterBase.h.

◆ _whatToFit

std::string lsst::jointcal::FitterBase::_whatToFit
protected

Definition at line 154 of file FitterBase.h.


The documentation for this class was generated from the following files: