LSSTApplications
20.0.0
LSSTDataManagementBasePackage
|
Class that handles the photometric least squares problem.
More...
#include <PhotometryFit.h>
Class that handles the photometric least squares problem.
Definition at line 45 of file PhotometryFit.h.
◆ PhotometryFit() [1/3]
Construct a photometry fitter.
- Parameters
-
associations | The associations catalog to use in the fitter. |
photometryModel | The model to build the fitter for. |
Definition at line 53 of file PhotometryFit.h.
57 _fittingFluxes(
false),
58 _photometryModel(photometryModel),
◆ PhotometryFit() [2/3]
lsst::jointcal::PhotometryFit::PhotometryFit |
( |
PhotometryFit const & |
| ) |
|
|
delete |
No copy or move: there is only ever one fitter of a given type.
◆ PhotometryFit() [3/3]
lsst::jointcal::PhotometryFit::PhotometryFit |
( |
PhotometryFit && |
| ) |
|
|
delete |
◆ assignIndices()
void lsst::jointcal::PhotometryFit::assignIndices |
( |
std::string const & |
whatToFit | ) |
|
|
overridevirtual |
Set parameters to fit and assign indices in the big matrix.
- Parameters
-
[in] | whatToFit | Valid strings : "Model", "Fluxes", which define which parameter sets are going to be fitted. whatToFit="Model Fluxes" will set both parameter sets variable when computing derivatives. Provided it contains "Model", whatToFit is passed over to the PhotometryModel, and can hence be used to control more finely which subsets of the photometric model are being fitted, if the the actual PhotometryModel implements such a possibility. |
Implements lsst::jointcal::FitterBase.
Definition at line 186 of file PhotometryFit.cc.
190 _fittingFluxes = (
_whatToFit.
find(
"Fluxes") != std::string::npos);
193 _nParModel = (_fittingModel) ? _photometryModel->assignIndices(whatToFit, 0) : 0;
196 if (_fittingFluxes) {
202 fittedStar->setIndexInMatrix(ipar);
206 _nParFluxes = ipar - _nParModel;
209 "nParameters total: " <<
_nParTot <<
" model: " << _nParModel <<
" fluxes: " << _nParFluxes);
◆ 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 ");
◆ getModel()
Return the model being fit.
Definition at line 88 of file PhotometryFit.h.
88 {
return _photometryModel; }
◆ leastSquareDerivatives()
void lsst::jointcal::FitterBase::leastSquareDerivatives |
( |
TripletList & |
tripletList, |
|
|
Eigen::VectorXd & |
grad |
|
) |
| const |
|
inherited |
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) {
◆ 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 = "" |
|
) |
| |
|
inherited |
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()
void lsst::jointcal::PhotometryFit::offsetParams |
( |
Eigen::VectorXd const & |
delta | ) |
|
|
overridevirtual |
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 |
Implements lsst::jointcal::FitterBase.
Definition at line 212 of file PhotometryFit.cc.
214 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
215 "PhotometryFit::offsetParams : the provided vector length is not compatible with "
216 "the current whatToFit setting");
217 if (_fittingModel) _photometryModel->offsetParams(delta);
219 if (_fittingFluxes) {
224 Eigen::Index index = fittedStar->getIndexInMatrix();
225 _photometryModel->offsetFittedStar(*fittedStar, delta(index));
◆ 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 | ) |
|
|
protectedinherited |
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 | ) |
|
|
protectedinherited |
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 |
|
virtualinherited |
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()
void lsst::jointcal::PhotometryFit::saveChi2MeasContributions |
( |
std::string const & |
filename | ) |
const |
|
overrideprotectedvirtual |
Save a CSV file containing residuals of measurement terms.
Implements lsst::jointcal::FitterBase.
Definition at line 230 of file PhotometryFit.cc.
234 ofile <<
"#id" << separator <<
"xccd" << separator <<
"yccd" << separator;
235 ofile <<
"mag" << separator <<
"instMag" << separator <<
"instMagErr" << separator;
236 ofile <<
"instFlux" << separator <<
"instFluxErr" << separator;
237 ofile <<
"inputFlux" << separator <<
"inputFluxErr" << separator;
238 ofile <<
"transformedFlux" << separator <<
"transformedFluxErr" << separator;
239 ofile <<
"fittedFlux" << separator;
240 ofile <<
"mjd" << separator <<
"color" << separator;
241 ofile <<
"fsindex" << separator;
242 ofile <<
"ra" << separator <<
"dec" << separator;
243 ofile <<
"chi2" << separator <<
"nm" << separator;
244 ofile <<
"chip" << separator <<
"visit" << separator <<
std::endl;
246 ofile <<
"#id in source catalog" << separator <<
"coordinates in CCD" << separator << separator;
247 ofile <<
"fitted magnitude" << separator <<
"measured magnitude" << separator
248 <<
"measured magnitude error" << separator;
249 ofile <<
"measured instrumental flux (ADU)" << separator <<
"measured instrument flux error" << separator;
250 ofile <<
"measured flux (nJy)" << separator <<
"measured flux error" << separator;
251 ofile << separator << separator;
252 ofile <<
"fitted flux (nJy)" << separator;
253 ofile <<
"modified Julian date of the measurement" << separator <<
"currently unused" << separator;
254 ofile <<
"unique index of the fittedStar" << separator;
255 ofile <<
"on-sky position of fitted star" << separator << separator;
256 ofile <<
"contribution to Chi2 (1 dof)" << separator <<
"number of measurements of this FittedStar"
258 ofile <<
"chip id" << separator <<
"visit id" <<
std::endl;
261 for (
auto const &ccdImage : ccdImageList) {
262 const MeasuredStarList &cat = ccdImage->getCatalogForFit();
263 for (
auto const &measuredStar : cat) {
264 if (!measuredStar->isValid())
continue;
266 double instFluxErr = _photometryModel->tweakFluxError(*measuredStar);
267 double flux = _photometryModel->transform(*ccdImage, *measuredStar);
268 double fluxErr = _photometryModel->transformError(*ccdImage, *measuredStar);
269 double jd = ccdImage->getMjd();
271 double residual = _photometryModel->computeResidual(*ccdImage, *measuredStar);
272 double chi2Val =
std::pow(residual / fluxErr, 2);
275 ofile << measuredStar->getId() << separator << measuredStar->x << separator << measuredStar->y
277 ofile << fittedStar->getMag() << separator << measuredStar->getInstMag() << separator
278 << measuredStar->getInstMagErr() << separator;
279 ofile << measuredStar->getInstFlux() << separator << instFluxErr << separator;
280 ofile << measuredStar->getFlux() << separator << measuredStar->getFluxErr() << separator;
281 ofile << flux << separator << fluxErr << separator << fittedStar->getFlux() << separator;
282 ofile << jd << separator << fittedStar->color << separator;
283 ofile << fittedStar->getIndexInMatrix() << separator;
284 ofile << fittedStar->x << separator << fittedStar->y << separator;
285 ofile << chi2Val << separator << fittedStar->getMeasurementCount() << separator;
286 ofile << ccdImage->getCcdId() << separator << ccdImage->getVisit() <<
std::endl;
◆ saveChi2RefContributions()
void lsst::jointcal::PhotometryFit::saveChi2RefContributions |
( |
std::string const & |
filename | ) |
const |
|
overrideprotectedvirtual |
Save a CSV file containing residuals of reference terms.
Implements lsst::jointcal::FitterBase.
Definition at line 291 of file PhotometryFit.cc.
295 ofile <<
"#ra" << separator <<
"dec " << separator;
296 ofile <<
"mag" << separator <<
"color" << separator;
297 ofile <<
"refFlux" << separator <<
"refFluxErr" << separator;
298 ofile <<
"fittedFlux" << separator <<
"fittedFluxErr" << separator;
299 ofile <<
"fsindex" << separator <<
"chi2" << separator <<
"nm" <<
std::endl;
301 ofile <<
"#coordinates of fittedStar" << separator << separator;
302 ofile <<
"magnitude" << separator <<
"currently unused" << separator;
303 ofile <<
"refStar flux (nJy)" << separator <<
"refStar fluxErr" << separator;
304 ofile <<
"fittedStar flux (nJy)" << separator <<
"fittedStar fluxErr" << separator;
305 ofile <<
"unique index of the fittedStar" << separator <<
"refStar contribution to Chi2 (1 dof)"
306 << separator <<
"number of measurements of this FittedStar" <<
std::endl;
309 const FittedStarList &fittedStarList =
_associations->fittedStarList;
310 for (
auto const &fittedStar : fittedStarList) {
311 const RefStar *refStar = fittedStar->getRefStar();
312 if (refStar ==
nullptr)
continue;
314 double chi2 =
std::pow(((fittedStar->getFlux() - refStar->getFlux()) / refStar->getFluxErr()), 2);
317 ofile << fittedStar->x << separator << fittedStar->y << separator;
318 ofile << fittedStar->getMag() << separator << fittedStar->color << separator;
319 ofile << refStar->getFlux() << separator << refStar->getFluxErr() << separator;
320 ofile << fittedStar->getFlux() << separator << fittedStar->getFluxErr() << separator;
321 ofile << fittedStar->getIndexInMatrix() << separator << chi2 << separator
322 << fittedStar->getMeasurementCount() <<
std::endl;
◆ _associations
◆ _lastNTrip
Eigen::Index lsst::jointcal::FitterBase::_lastNTrip |
|
protectedinherited |
◆ _log
◆ _nMeasuredStars
Eigen::Index lsst::jointcal::FitterBase::_nMeasuredStars |
|
protectedinherited |
◆ _nParTot
Eigen::Index lsst::jointcal::FitterBase::_nParTot |
|
protectedinherited |
◆ _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/PhotometryFit.h
- /j/snowflake/release/lsstsw/stack/1a1d771/Linux64/jointcal/20.0.0/src/PhotometryFit.cc
T setprecision(T... args)
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.
std::list< std::shared_ptr< CcdImage > > CcdImageList
#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.
FitterBase(std::shared_ptr< Associations > associations)
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()
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
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.