LSSTApplications
20.0.0
LSSTDataManagementBasePackage
|
Base class for fitters.
More...
#include <FitterBase.h>
|
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. 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...
|
|
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.
◆ FitterBase() [1/3]
◆ 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 |
◆ accumulateStatImageList()
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 |
◆ computeChi2()
Returns the chi2 for the current state.
Definition at line 42 of file FitterBase.cc.
◆ findOutliers()
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] | nSigmaCut | Number of sigma to select on. |
[out] | msOutliers | list of MeasuredStar outliers to populate |
[out] | fsOutliers | list of FittedStar outliers to populate |
- Returns
- Total number of outliers that were removed.
Definition at line 52 of file FitterBase.cc.
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;
77 Eigen::VectorXi affectedParams(
_nParTot);
78 affectedParams.setZero();
82 for (
auto chi2 = chi2List.rbegin(); chi2 != chi2List.rend(); ++chi2) {
83 if (chi2->chi2 < cut)
break;
88 auto measuredStar = std::dynamic_pointer_cast<MeasuredStar>(chi2->star);
90 if (measuredStar ==
nullptr) {
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);
99 indices.
push_back(fittedStar->getIndexInMatrix());
100 LOGLS_TRACE(
_log,
"Removing refStar " << *(fittedStar->getRefStar()) <<
" chi2: " << chi2->chi2);
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: "
112 LOGLS_TRACE(
_log,
"Removing measStar " << *measuredStar <<
" chi2: " << chi2->chi2);
120 for (
auto const &i : indices) {
121 if (affectedParams(i) != 0) {
128 if (measuredStar ==
nullptr) {
130 fsOutliers.push_back(fittedStar);
133 msOutliers.push_back(measuredStar);
136 for (
auto const &i : indices) {
142 LOGLS_INFO(
_log,
"findOutliers: found " << msOutliers.size() <<
" meas outliers and " << fsOutliers.size()
143 <<
" ref outliers ");
◆ 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
-
tripletList | tripletList of (row,col,value) representing the Jacobian of the chi2. |
grad | The gradient of the chi2. |
Definition at line 317 of file FitterBase.cc.
319 for (
auto const &ccdImage : ccdImageList) {
◆ 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] | whatToFit | See child method assignIndices for valid string values. |
[in] | nSigmaCut | How many sigma to reject outliers at. Outlier rejection ignored for nSigmaCut=0. |
[in] | doRankUpdate | Use 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] | doLineSearch | Use 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] | dumpMatrixFile | Write 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.
179 TripletList tripletList(nTrip);
188 LOGLS_DEBUG(
_log,
"End of triplet filling, ntrip = " << tripletList.size());
194 << hessian.rows() <<
" non-zeros=" << hessian.nonZeros()
195 <<
" filling-frac = " << hessian.nonZeros() /
std::pow(hessian.rows(), 2));
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());
202 dumpMatrixAndGradient(hessian, grad, dumpMatrixFile,
_log);
207 if (chol.info() != Eigen::Success) {
217 Eigen::VectorXd delta = chol.solve(grad);
219 scale = _lineSearch(delta);
225 LOGL_ERROR(
_log,
"chi2 is not finite. Aborting outlier rejection.");
229 if (currentChi2.chi2 > oldChi2 && totalMeasOutliers + totalRefOutliers != 0) {
230 LOGL_WARN(
_log,
"chi2 went up, skipping outlier rejection loop");
234 oldChi2 = currentChi2.chi2;
236 if (nSigmaCut == 0)
break;
237 MeasuredStarList msOutliers;
238 FittedStarList fsOutliers;
241 totalMeasOutliers += msOutliers.size();
242 totalRefOutliers += fsOutliers.size();
243 if (nOutliers == 0)
break;
244 TripletList outlierTriplets(nOutliers);
254 H.setFromTriplets(outlierTriplets.begin(), outlierTriplets.end());
255 chol.update(H,
false );
266 LOGLS_DEBUG(
_log,
"Triplets recomputed, ntrip = " << nextTripletList.size());
268 hessian = createHessian(
_nParTot, nextTripletList);
269 nextTripletList.clear();
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) {
284 if (nSigmaCut != 0) {
285 LOGLS_INFO(
_log,
"Number of outliers (Measured + Reference = Total): "
286 << totalMeasOutliers <<
" + " << totalRefOutliers <<
" = "
287 << totalMeasOutliers + totalRefOutliers);
◆ 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] | delta | vector of offsets to apply |
Implemented in lsst::jointcal::AstrometryFit, and lsst::jointcal::PhotometryFit.
◆ operator=() [1/2]
◆ operator=() [2/2]
◆ outliersContributions()
Contributions to derivatives from (presumably) outlier terms.
No discarding done.
Definition at line 292 of file FitterBase.cc.
294 for (
auto &outlier : msOutliers) {
295 MeasuredStarList tmp;
296 tmp.push_back(outlier);
297 const CcdImage &ccdImage = outlier->getCcdImage();
◆ 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.
304 for (
auto &measuredStar : outliers) {
305 auto fittedStar = measuredStar->getFittedStar();
306 measuredStar->setValid(
false);
307 fittedStar->getMeasurementCount()--;
◆ 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.
312 for (
auto &fittedStar : outliers) {
313 fittedStar->setRefStar(
nullptr);
◆ 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.
327 auto pos = baseName.
find(replaceStr);
329 measFilename.replace(pos, replaceStr.
size(),
"-meas.csv");
331 refFilename.replace(pos, replaceStr.
size(),
"-ref.csv");
◆ saveChi2MeasContributions()
virtual void lsst::jointcal::FitterBase::saveChi2MeasContributions |
( |
std::string const & |
filename | ) |
const |
|
protectedpure virtual |
◆ saveChi2RefContributions()
virtual void lsst::jointcal::FitterBase::saveChi2RefContributions |
( |
std::string const & |
filename | ) |
const |
|
protectedpure virtual |
◆ _associations
◆ _lastNTrip
Eigen::Index lsst::jointcal::FitterBase::_lastNTrip |
|
protected |
◆ _log
◆ _nMeasuredStars
Eigen::Index lsst::jointcal::FitterBase::_nMeasuredStars |
|
protected |
◆ _nParTot
Eigen::Index lsst::jointcal::FitterBase::_nParTot |
|
protected |
◆ _whatToFit
The documentation for this class was generated from the following files:
- /j/snowflake/release/lsstsw/stack/1a1d771/Linux64/jointcal/20.0.0/include/lsst/jointcal/FitterBase.h
- /j/snowflake/release/lsstsw/stack/1a1d771/Linux64/jointcal/20.0.0/src/FitterBase.cc
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.
#define LOGL_WARN(logger, message...)
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.
#define LOGLS_INFO(logger, message)
#define LOGLS_ERROR(logger, message)
#define LOGLS_TRACE(logger, message)
virtual void accumulateStatRefStars(Chi2Accumulator &accum) const =0
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars.
#define LOGLS_WARN(logger, message)
#define LOGL_ERROR(logger, message...)
std::size_t findOutliers(double nSigmaCut, MeasuredStarList &msOutliers, FittedStarList &fsOutliers) const
Find Measurements and references contributing more than a cut, computed as.
std::shared_ptr< Associations > _associations
void leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const
Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting.
virtual void getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, IndexVector &indices) const =0
Set the indices of a measured star from the full matrix, for outlier removal.
MinimizeResult
Return value of minimize()
virtual void saveChi2MeasContributions(std::string const &filename) const =0
Save a CSV file containing residuals of measurement terms.
#define LOGLS_DEBUG(logger, message)
void outliersContributions(MeasuredStarList &msOutliers, FittedStarList &fsOutliers, TripletList &tripletList, Eigen::VectorXd &grad)
Contributions to derivatives from (presumably) outlier terms.
Chi2Statistic computeChi2() const
Returns the chi2 for the current state.
virtual void assignIndices(std::string const &whatToFit)=0
Set parameters to fit and assign indices in the big matrix.
void removeMeasOutliers(MeasuredStarList &outliers)
Remove measuredStar outliers from the fit. No Refit done.
Eigen::SparseMatrix< double, 0, Eigen::Index > SparseMatrixD
virtual void offsetParams(Eigen::VectorXd const &delta)=0
Offset the parameters by the requested quantities.
def scale(algorithm, min, max=None, frame=None)
void removeRefOutliers(FittedStarList &outliers)
Remove refStar outliers from the fit. No Refit done.
Eigen::Index _nMeasuredStars
virtual void saveChi2RefContributions(std::string const &filename) const =0
Save a CSV file containing residuals of reference terms.