LSST Applications g0f08755f38+9c285cab97,g1635faa6d4+bcae251498,g1653933729+a8ce1bb630,g1a0ca8cf93+bf6eb00ceb,g28da252d5a+0829b12dee,g29321ee8c0+18ecbd06b3,g2bbee38e9b+9634bc57db,g2bc492864f+9634bc57db,g2cdde0e794+c2c89b37c4,g3156d2b45e+41e33cbcdc,g347aa1857d+9634bc57db,g35bb328faa+a8ce1bb630,g3a166c0a6a+9634bc57db,g3e281a1b8c+9f2c4e2fc3,g414038480c+077ccc18e7,g41af890bb2+fde0dd39b6,g5fbc88fb19+17cd334064,g7642f7d749+9c285cab97,g781aacb6e4+a8ce1bb630,g80478fca09+55a9465950,g82479be7b0+ed77629bff,g858d7b2824+9c285cab97,g9125e01d80+a8ce1bb630,g9726552aa6+10f999ec6a,ga5288a1d22+2a84bb7594,gacf8899fa4+c69c5206e8,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gbd46683f8f+1c79523530,gc28159a63d+9634bc57db,gcf0d15dbbd+4b7d09cae4,gda3e153d99+9c285cab97,gda6a2b7d83+4b7d09cae4,gdaeeff99f8+1711a396fd,ge2409df99d+dfd3d5294a,ge79ae78c31+9634bc57db,gf0baf85859+147a0692ba,gf3967379c6+02b11634a5,w.2024.46
LSST Data Management Base Package
|
Class that handles the photometric least squares problem. More...
#include <PhotometryFit.h>
Public Member Functions | |
PhotometryFit (std::shared_ptr< Associations > associations, std::shared_ptr< PhotometryModel > photometryModel) | |
Construct a photometry fitter. | |
PhotometryFit (PhotometryFit const &)=delete | |
No copy or move: there is only ever one fitter of a given type. | |
PhotometryFit (PhotometryFit &&)=delete | |
PhotometryFit & | operator= (PhotometryFit const &)=delete |
PhotometryFit & | operator= (PhotometryFit &&)=delete |
void | assignIndices (std::string const &whatToFit) override |
Set parameters to fit and assign indices in the big matrix. | |
void | offsetParams (Eigen::VectorXd const &delta) override |
Offset the parameters by the requested quantities. | |
std::shared_ptr< PhotometryModel > | getModel () const |
Return the model being fit. | |
MinimizeResult | minimize (std::string const &whatToFit, double nSigmaCut=0, double sigmaRelativeTolerance=0, bool doRankUpdate=true, bool doLineSearch=false, std::string const &dumpMatrixFile="") |
Does a 1 step minimization, assuming a linear model. | |
Chi2Statistic | computeChi2 () const |
Returns the chi2 for the current state. | |
void | leastSquareDerivatives (TripletList &tripletList, Eigen::VectorXd &grad) const |
Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting. | |
virtual void | saveChi2Contributions (std::string const &baseName) const |
Save the full chi2 term per star that was used in the minimization, for debugging. | |
Protected Member Functions | |
void | saveChi2MeasContributions (std::string const &filename) const override |
Save a CSV file containing residuals of measurement terms. | |
void | saveChi2RefContributions (std::string const &filename) const override |
Save a CSV file containing residuals of reference terms. | |
std::size_t | findOutliers (double nSigmaCut, MeasuredStarList &msOutliers, FittedStarList &fsOutliers, double &cut) const |
Find Measurements and references contributing more than a cut, computed as. | |
void | outliersContributions (MeasuredStarList &msOutliers, FittedStarList &fsOutliers, TripletList &tripletList, Eigen::VectorXd &grad) |
Contributions to derivatives from (presumably) outlier terms. | |
void | removeMeasOutliers (MeasuredStarList &outliers) |
Remove measuredStar outliers from the fit. No Refit done. | |
void | removeRefOutliers (FittedStarList &outliers) |
Remove refStar outliers from the fit. No Refit done. | |
Protected Attributes | |
std::shared_ptr< Associations > | _associations |
std::string | _whatToFit |
Eigen::Index | _lastNTrip |
Eigen::Index | _nTotal |
Eigen::Index | _nModelParams |
Eigen::Index | _nStarParams |
LOG_LOGGER | _log |
Private Member Functions | |
void | accumulateStatImageList (CcdImageList const &ccdImageList, Chi2Accumulator &accum) const override |
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for measurements. | |
void | accumulateStatRefStars (Chi2Accumulator &accum) const override |
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars. | |
void | getIndicesOfMeasuredStar (MeasuredStar const &measuredStar, IndexVector &indices) const override |
this routine is to be used only in the framework of outlier removal | |
void | leastSquareDerivativesMeasurement (CcdImage const &ccdImage, TripletList &tripletList, Eigen::VectorXd &grad, MeasuredStarList const *measuredStarList=nullptr) const override |
Compute the derivatives of the measured stars and model for one CcdImage. | |
void | leastSquareDerivativesReference (FittedStarList const &fittedStarList, TripletList &tripletList, Eigen::VectorXd &grad) const override |
Compute the derivatives of the reference terms. | |
Class that handles the photometric least squares problem.
Definition at line 46 of file PhotometryFit.h.
|
inline |
Construct a photometry fitter.
associations | The associations catalog to use in the fitter. |
photometryModel | The model to build the fitter for. |
Definition at line 54 of file PhotometryFit.h.
|
delete |
No copy or move: there is only ever one fitter of a given type.
|
delete |
|
overrideprivatevirtual |
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for measurements.
Implements lsst::jointcal::FitterBase.
Definition at line 134 of file PhotometryFit.cc.
|
overrideprivatevirtual |
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars.
Implements lsst::jointcal::FitterBase.
Definition at line 153 of file PhotometryFit.cc.
|
overridevirtual |
Set parameters to fit and assign indices in the big matrix.
[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 185 of file PhotometryFit.cc.
|
inherited |
Returns the chi2 for the current state.
Definition at line 43 of file FitterBase.cc.
|
protectedinherited |
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.
[in] | nSigmaCut | Number of sigma to select on. |
[out] | msOutliers | list of MeasuredStar outliers to populate |
[out] | fsOutliers | list of FittedStar outliers to populate |
[out] | cut | value of chi2 that defines which objects are outliers |
Definition at line 53 of file FitterBase.cc.
|
overrideprivatevirtual |
this routine is to be used only in the framework of outlier removal
it fills the array of indices of parameters that a Measured star constrains. Not really all of them if you check.
Implements lsst::jointcal::FitterBase.
Definition at line 173 of file PhotometryFit.cc.
|
inline |
Return the model being fit.
Definition at line 87 of file PhotometryFit.h.
|
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.
tripletList | tripletList of (row,col,value) representing the Jacobian of the chi2. |
grad | The gradient of the chi2. |
Definition at line 347 of file FitterBase.cc.
|
overrideprivatevirtual |
Compute the derivatives of the measured stars and model for one CcdImage.
The last argument will process a sub-list for outlier removal.
Implements lsst::jointcal::FitterBase.
Definition at line 45 of file PhotometryFit.cc.
|
overrideprivatevirtual |
Compute the derivatives of the reference terms.
Implements lsst::jointcal::FitterBase.
Definition at line 97 of file PhotometryFit.cc.
|
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.
[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] | sigmaRelativeTolerance | Percentage change in the chi2 cut for outliers tolerated for termination. If value is zero, minimization iterations will continue until there are no outliers. |
[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)
|
Definition at line 179 of file FitterBase.cc.
|
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.
[in] | delta | vector of offsets to apply |
Implements lsst::jointcal::FitterBase.
Definition at line 211 of file PhotometryFit.cc.
|
delete |
|
delete |
|
protectedinherited |
Contributions to derivatives from (presumably) outlier terms.
No discarding done.
Definition at line 322 of file FitterBase.cc.
|
protectedinherited |
Remove measuredStar outliers from the fit. No Refit done.
Definition at line 333 of file FitterBase.cc.
|
protectedinherited |
Remove refStar outliers from the fit. No Refit done.
Definition at line 341 of file FitterBase.cc.
|
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 355 of file FitterBase.cc.
|
overrideprotectedvirtual |
Save a CSV file containing residuals of measurement terms.
Implements lsst::jointcal::FitterBase.
Definition at line 229 of file PhotometryFit.cc.
|
overrideprotectedvirtual |
Save a CSV file containing residuals of reference terms.
Implements lsst::jointcal::FitterBase.
Definition at line 289 of file PhotometryFit.cc.
|
protectedinherited |
Definition at line 165 of file FitterBase.h.
|
protectedinherited |
Definition at line 168 of file FitterBase.h.
|
protectedinherited |
Definition at line 174 of file FitterBase.h.
|
protectedinherited |
Definition at line 170 of file FitterBase.h.
|
protectedinherited |
Definition at line 171 of file FitterBase.h.
|
protectedinherited |
Definition at line 169 of file FitterBase.h.
|
protectedinherited |
Definition at line 166 of file FitterBase.h.