LSST Applications g063fba187b+cac8b7c890,g0f08755f38+6aee506743,g1653933729+a8ce1bb630,g168dd56ebc+a8ce1bb630,g1a2382251a+b4475c5878,g1dcb35cd9c+8f9bc1652e,g20f6ffc8e0+6aee506743,g217e2c1bcf+73dee94bd0,g28da252d5a+1f19c529b9,g2bbee38e9b+3f2625acfc,g2bc492864f+3f2625acfc,g3156d2b45e+6e55a43351,g32e5bea42b+1bb94961c2,g347aa1857d+3f2625acfc,g35bb328faa+a8ce1bb630,g3a166c0a6a+3f2625acfc,g3e281a1b8c+c5dd892a6c,g3e8969e208+a8ce1bb630,g414038480c+5927e1bc1e,g41af890bb2+8a9e676b2a,g7af13505b9+809c143d88,g80478fca09+6ef8b1810f,g82479be7b0+f568feb641,g858d7b2824+6aee506743,g89c8672015+f4add4ffd5,g9125e01d80+a8ce1bb630,ga5288a1d22+2903d499ea,gb58c049af0+d64f4d3760,gc28159a63d+3f2625acfc,gcab2d0539d+b12535109e,gcf0d15dbbd+46a3f46ba9,gda6a2b7d83+46a3f46ba9,gdaeeff99f8+1711a396fd,ge79ae78c31+3f2625acfc,gef2f8181fd+0a71e47438,gf0baf85859+c1f95f4921,gfa517265be+6aee506743,gfa999e8aa5+17cd334064,w.2024.51
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
lsst::jointcal::AstrometryFit Class Reference

Class that handles the astrometric least squares problem. More...

#include <AstrometryFit.h>

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

Public Member Functions

 AstrometryFit (std::shared_ptr< Associations > associations, std::shared_ptr< AstrometryModel > astrometryModel, double posError)
 this is the only constructor
 
 AstrometryFit (AstrometryFit const &)=delete
 No copy or move: there is only ever one fitter of a given type.
 
 AstrometryFit (AstrometryFit &&)=delete
 
AstrometryFitoperator= (AstrometryFit const &)=delete
 
AstrometryFitoperator= (AstrometryFit &&)=delete
 
void assignIndices (std::string const &whatToFit) override
 Set parameters to fit and assign indices in the big matrix.
 
void freezeErrorTransform ()
 The transformations used to propagate errors are freezed to the current state.
 
void offsetParams (Eigen::VectorXd const &delta) override
 Offset the parameters by the requested quantities.
 
std::shared_ptr< AstrometryModelgetModel () const
 Return the model being fit.
 
void checkStuff ()
 DEBUGGING routine.
 
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 leastSquareDerivativesMeasurement (CcdImage const &ccdImage, TripletList &tripletList, Eigen::VectorXd &grad, MeasuredStarList const *msList=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.
 
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
 

Detailed Description

Class that handles the astrometric least squares problem.

This is the class that actually computes the quantities required to carry out a LS astrometric fit wrt distortion mappings and coordinates of common objects. Namely it computes the Jacobian and gradient of the chi2 (w.r.t. parameters), and the Chi2 itself. It interfaces with the actual modelling of distortions via a mimimum virtual interface AstrometryModel, and the actual mappings via an other virtual interface : Mapping.

In short AstrometryFit aims at computing derivatives of least quares. The terms of the chi2 are of two kinds:

kind 1 -> (T(X_M) - p(F))^T W (T(X_M) - p(F))

with X_M is a measured (2d) position in CCD coordinates, F refers to the position of the object in some space, defined in practise by p. There is one such term per measurement. The default setup would be that p is the projection from sky to some tangent plane and hence T maps the CCD coordinates onto this TP. p is obtained via the DistorsionModel and can be different for all CcdImage's. Depending on what is beeing fitted, one could imagine cases where the projector p is the same for all CcdImages.

Kind 2 -> (p'(F)-p'(R))^T W_R (p'(F)-p'(R)) R refers to some externally-provided reference object position, and p' to some projector from sky to some plane. The reference objects define the overall coordinate frame, which is required when all T and all F are fitted simultaneously. There is one such term per external reference object. There can be more F (fitted) objects than R (reference) objects.

In the same framework, one can fit relative transforms between images by setting p = Identity for all input CcdImages and not fitting T for one of the CcdImage's. One does not need reference object and would then naturally not have any Kind 2 terms.

Definition at line 78 of file AstrometryFit.h.

Constructor & Destructor Documentation

◆ AstrometryFit() [1/3]

lsst::jointcal::AstrometryFit::AstrometryFit ( std::shared_ptr< Associations > associations,
std::shared_ptr< AstrometryModel > astrometryModel,
double posError )

this is the only constructor

Definition at line 66 of file AstrometryFit.cc.

68 : FitterBase(std::move(associations)),
69 _astrometryModel(std::move(astrometryModel)),
70 _epoch(_associations->getEpoch()),
71 _posError(posError) {
72 _log = LOG_GET("lsst.jointcal.AstrometryFit");
73}
#define LOG_GET(logger)
Returns a Log object associated with logger.
Definition Log.h:75
std::shared_ptr< Associations > _associations
Definition FitterBase.h:165
FitterBase(std::shared_ptr< Associations > associations)
Definition FitterBase.h:57
T move(T... args)

◆ AstrometryFit() [2/3]

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

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

◆ AstrometryFit() [3/3]

lsst::jointcal::AstrometryFit::AstrometryFit ( AstrometryFit && )
delete

Member Function Documentation

◆ accumulateStatImageList()

void lsst::jointcal::AstrometryFit::accumulateStatImageList ( CcdImageList const & ccdImageList,
Chi2Accumulator & accum ) const
overrideprivatevirtual

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

Implements lsst::jointcal::FitterBase.

Definition at line 372 of file AstrometryFit.cc.

372 {
373 for (auto const &ccdImage : ccdImageList) {
374 accumulateStatImage(*ccdImage, accum);
375 }
376}

◆ accumulateStatRefStars()

void lsst::jointcal::AstrometryFit::accumulateStatRefStars ( Chi2Accumulator & accum) const
overrideprivatevirtual

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

Note
the math in this method and leastSquareDerivativesReference() must be kept consistent, in terms of +/- convention, definition of model, etc.

Implements lsst::jointcal::FitterBase.

Definition at line 378 of file AstrometryFit.cc.

378 {
379 /**********************************************************************/
382 /**********************************************************************/
383
384 /* If you wonder why we project here, read comments in
385 AstrometryFit::leastSquareDerivativesReference(TripletList &TList, Eigen::VectorXd &Rhs) */
386 FittedStarList &fittedStarList = _associations->fittedStarList;
387 TanRaDecToPixel proj(AstrometryTransformLinear(), Point(0., 0.));
388 for (auto const &fs : fittedStarList) {
389 auto rs = fs->getRefStar();
390 if (rs == nullptr) continue;
391 proj.setTangentPoint(*fs);
392 // fs projects to (0,0), no need to compute its transform.
393 FatPoint rsProj;
394 proj.transformPosAndErrors(*rs, rsProj);
395 // TO DO : account for proper motions.
396 double chi2 = computeProjectedRefStarChi2(rsProj);
397 accum.addEntry(chi2, 2, fs);
398 }
399}

◆ assignIndices()

void lsst::jointcal::AstrometryFit::assignIndices ( std::string const & whatToFit)
overridevirtual

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

Parameters
[in]whatToFitValid strings: zero or more of "Distortions", "Positions", "PM" which define which parameter set is going to be variable when computing derivatives (leastSquareDerivatives) and minimizing (minimize()). whatToFit="Positions Distortions" will minimize w.r.t mappings and objects positions, and not w.r.t proper motions and refraction modeling. However if proper motions and/or refraction parameters have already been set, then they are accounted for when computing residuals. The string is forwarded to the AstrometryModel, and it can then be used to turn subsets of distortion parameter on or off, if the AstrometryModel implements such a thing.

Implements lsst::jointcal::FitterBase.

Definition at line 423 of file AstrometryFit.cc.

423 {
424 _whatToFit = whatToFit;
425 LOGLS_INFO(_log, "assignIndices: Now fitting " << whatToFit);
426 _fittingDistortions = (_whatToFit.find("Distortions") != std::string::npos);
427 _fittingPos = (_whatToFit.find("Positions") != std::string::npos);
428 _fittingPM = (_whatToFit.find("PM") != std::string::npos);
429 // When entering here, we assume that whatToFit has already been interpreted.
430
431 _nModelParams = 0;
432 if (_fittingDistortions) _nModelParams = _astrometryModel->assignIndices(_whatToFit, 0);
434
435 if (_fittingPos) {
436 FittedStarList &fittedStarList = _associations->fittedStarList;
437 for (auto &fittedStar : fittedStarList) {
438 // the parameter layout here is used also
439 // - when filling the derivatives
440 // - when updating (offsetParams())
441 // - in GetMeasuredStarIndices
442 fittedStar->setIndexInMatrix(ipar);
443 ipar += 2;
444 // TODO: left as reference for when we implement PM fitting
445 // if ((_fittingPM)&fittedStar->mightMove) ipar += 2;
446 }
447 }
449 _nTotal = ipar;
450 LOGLS_DEBUG(_log, "nParameters total: " << _nTotal << " model: " << _nModelParams
451 << " values: " << _nStarParams);
452}
#define LOGLS_INFO(logger, message)
Log a info-level message using an iostream-based interface.
Definition Log.h:639
#define LOGLS_DEBUG(logger, message)
Log a debug-level message using an iostream-based interface.
Definition Log.h:619
T find(T... args)

◆ checkStuff()

void lsst::jointcal::AstrometryFit::checkStuff ( )

DEBUGGING routine.

Definition at line 480 of file AstrometryFit.cc.

480 {
481 const char *what2fit[] = {"Positions", "Distortions", "Positions Distortions"};
482 // DEBUG
483 for (unsigned k = 0; k < sizeof(what2fit) / sizeof(what2fit[0]); ++k) {
484 assignIndices(what2fit[k]);
485 TripletList tripletList(10000);
486 Eigen::VectorXd grad(_nTotal);
487 grad.setZero();
488 leastSquareDerivatives(tripletList, grad);
489 SparseMatrixD jacobian(_nTotal, tripletList.getNextFreeIndex());
490 jacobian.setFromTriplets(tripletList.begin(), tripletList.end());
491 SparseMatrixD hessian = jacobian * jacobian.transpose();
492 LOGLS_DEBUG(_log, "npar : " << _nTotal << ' ' << _nModelParams);
493 }
494}
Eigen::SparseMatrix< double, 0, Eigen::Index > SparseMatrixD
Definition Eigenstuff.h:35
void assignIndices(std::string const &whatToFit) override
Set parameters to fit and assign indices in the big matrix.
void leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const
Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting.

◆ computeChi2()

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

Returns the chi2 for the current state.

Definition at line 43 of file FitterBase.cc.

43 {
44 Chi2Statistic chi2;
45 accumulateStatImageList(_associations->getCcdImageList(), chi2);
47 // chi2.ndof contains the number of squares.
48 // So subtract the number of parameters.
49 chi2.ndof -= _nTotal;
50 return chi2;
51}
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.

◆ findOutliers()

std::size_t lsst::jointcal::FitterBase::findOutliers ( double nSigmaCut,
MeasuredStarList & msOutliers,
FittedStarList & fsOutliers,
double & cut ) const
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.

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

Definition at line 53 of file FitterBase.cc.

54 {
55 // collect chi2 contributions
56 Chi2List chi2List;
57 chi2List.reserve(_associations->getMaxMeasuredStars() + _associations->refStarList.size());
58 // contributions from measurement terms:
59 accumulateStatImageList(_associations->ccdImageList, chi2List);
60 // and from reference terms
61 accumulateStatRefStars(chi2List);
62
63 // compute some statistics
64 size_t nval = chi2List.size();
65 if (nval == 0) return 0;
66 sort(chi2List.begin(), chi2List.end());
67 double median = (nval & 1) ? chi2List[nval / 2].chi2
68 : 0.5 * (chi2List[nval / 2 - 1].chi2 + chi2List[nval / 2].chi2);
69 auto averageAndSigma = chi2List.computeAverageAndSigma();
70 LOGLS_DEBUG(_log, "findOutliers chi2 stat: mean/median/sigma " << averageAndSigma.first << '/' << median
71 << '/' << averageAndSigma.second);
72 cut = averageAndSigma.first + nSigmaCut * averageAndSigma.second;
73 /* For each of the parameters, we will not remove more than 1
74 measurement that contributes to constraining it. Keep track using
75 of what we are touching using an integer vector. This is the
76 trick that Marc Betoule came up to for outlier removals in "star
77 flats" fits. */
78 Eigen::VectorXi affectedParams(_nTotal);
79 affectedParams.setZero();
80
81 std::size_t nOutliers = 0; // returned to the caller
82 // start from the strongest outliers.
83 for (auto chi2 = chi2List.rbegin(); chi2 != chi2List.rend(); ++chi2) {
84 if (chi2->chi2 < cut) break; // because the array is sorted.
85 IndexVector indices;
86 /* now, we want to get the indices of the parameters this chi2
87 term depends on. We have to figure out which kind of term it
88 is; we use for that the type of the star attached to the Chi2Star. */
89 auto measuredStar = std::dynamic_pointer_cast<MeasuredStar>(chi2->star);
90 std::shared_ptr<FittedStar> fittedStar; // To add to fsOutliers if it is a reference outlier.
91 if (measuredStar == nullptr) {
92 // it is a reference outlier
93 fittedStar = std::dynamic_pointer_cast<FittedStar>(chi2->star);
94 if (fittedStar->getMeasurementCount() == 0) {
95 LOGLS_WARN(_log, "FittedStar with no measuredStars found as an outlier: "
96 << *fittedStar << " chi2: " << chi2->chi2);
97 continue;
98 }
99 if (_nStarParams == 0) {
101 "RefStar is outlier but not removed when not fitting FittedStar-RefStar values: "
102 << *(fittedStar->getRefStar()) << " chi2: " << chi2->chi2);
103 continue;
104 }
105 // NOTE: Stars contribute twice to astrometry (x,y), but once to photometry (flux),
106 // NOTE: but we only need to mark one index here because both will be removed with that star.
107 indices.push_back(fittedStar->getIndexInMatrix());
108 LOGLS_TRACE(_log, "Removing refStar " << *(fittedStar->getRefStar()) << " chi2: " << chi2->chi2);
109 /* One might think it would be useful to account for PM
110 parameters here, but it is just useless */
111 } else {
112 // it is a measurement outlier
113 auto tempFittedStar = measuredStar->getFittedStar();
114 if (tempFittedStar->getMeasurementCount() == 1 && tempFittedStar->getRefStar() == nullptr) {
115 LOGLS_WARN(_log, "FittedStar with 1 measuredStar and no refStar found as an outlier: "
116 << *tempFittedStar);
117 continue;
118 }
119 getIndicesOfMeasuredStar(*measuredStar, indices);
120 LOGLS_TRACE(_log, "Removing measStar " << *measuredStar << " chi2: " << chi2->chi2);
121 }
122
123 /* Find out if we already discarded a stronger outlier
124 constraining some parameter this one constrains as well. If
125 yes, we keep this one, because this stronger outlier could be
126 causing the large chi2 we have in hand. */
127 bool drop_it = true;
128 for (auto const &i : indices) {
129 if (affectedParams(i) != 0) {
130 drop_it = false;
131 }
132 }
133
134 if (drop_it) // store the outlier in one of the lists:
135 {
136 if (measuredStar == nullptr) {
137 // reference term
138 fsOutliers.push_back(fittedStar);
139 } else {
140 // measurement term
141 msOutliers.push_back(measuredStar);
142 }
143 // mark the parameters as directly changed when we discard this chi2 term.
144 for (auto const &i : indices) {
145 affectedParams(i)++;
146 }
147 nOutliers++;
148 }
149 } // end loop on measurements/references
150 LOGLS_DEBUG(_log, "findOutliers: found " << msOutliers.size() << " meas outliers and "
151 << fsOutliers.size() << " ref outliers ");
152
153 return nOutliers;
154}
#define LOGLS_WARN(logger, message)
Log a warn-level message using an iostream-based interface.
Definition Log.h:659
#define LOGLS_TRACE(logger, message)
Log a trace-level message using an iostream-based interface.
Definition Log.h:599
virtual void getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, IndexVector &indices) const =0
Set the indices of a measured star from the full matrix, for outlier removal.
T push_back(T... args)
T sort(T... args)

◆ freezeErrorTransform()

void lsst::jointcal::AstrometryFit::freezeErrorTransform ( )
inline

The transformations used to propagate errors are freezed to the current state.

The routine can be called when the mappings are roughly in place. After the call, the transformations used to propage errors are no longer affected when updating the mappings. This allows to have an exactly linear fit, which can be useful.

Definition at line 117 of file AstrometryFit.h.

117{ _astrometryModel->freezeErrorTransform(); }

◆ getIndicesOfMeasuredStar()

void lsst::jointcal::AstrometryFit::getIndicesOfMeasuredStar ( MeasuredStar const & measuredStar,
IndexVector & indices ) const
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 404 of file AstrometryFit.cc.

404 {
405 if (_fittingDistortions) {
406 const AstrometryMapping *mapping = _astrometryModel->getMapping(measuredStar.getCcdImage());
407 mapping->getMappingIndices(indices);
408 }
409 std::shared_ptr<FittedStar const> const fs = measuredStar.getFittedStar();
410 Eigen::Index fsIndex = fs->getIndexInMatrix();
411 if (_fittingPos) {
412 indices.push_back(fsIndex);
413 indices.push_back(fsIndex + 1);
414 }
415 // For securing the outlier removal, the next block is just useless
416 if (_fittingPM) {
417 for (std::size_t k = 0; k < 2; ++k) indices.push_back(fsIndex + 2 + k);
418 }
419 /* Should not put the index of refaction stuff or we will not be
420 able to remove more than 1 star at a time. */
421}

◆ getModel()

std::shared_ptr< AstrometryModel > lsst::jointcal::AstrometryFit::getModel ( ) const
inline

Return the model being fit.

Definition at line 122 of file AstrometryFit.h.

122{ return _astrometryModel; }

◆ 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
tripletListtripletList of (row,col,value) representing the Jacobian of the chi2.
gradThe gradient of the chi2.

Definition at line 347 of file FitterBase.cc.

347 {
348 auto ccdImageList = _associations->getCcdImageList();
349 for (auto const &ccdImage : ccdImageList) {
350 leastSquareDerivativesMeasurement(*ccdImage, tripletList, grad);
351 }
352 leastSquareDerivativesReference(_associations->fittedStarList, tripletList, grad);
353}
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.

◆ leastSquareDerivativesMeasurement()

void lsst::jointcal::AstrometryFit::leastSquareDerivativesMeasurement ( CcdImage const & ccdImage,
TripletList & tripletList,
Eigen::VectorXd & grad,
MeasuredStarList const * measuredStarList = nullptr ) const
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 107 of file AstrometryFit.cc.

109 {
110 /**********************************************************************/
111 /* @note the math in this method and accumulateStatImage() must be kept consistent,
112 * in terms of +/- convention, definition of model, etc. */
113 /**********************************************************************/
114
115 /* Setup */
116 /* this routine works in two different ways: either providing the
117 ccdImage, of providing the MeasuredStarList. In the latter case, the
118 ccdImage should match the one(s) in the list. */
119 if (msList) assert(&(msList->front()->getCcdImage()) == &ccdImage);
120
121 // get the Mapping
122 const AstrometryMapping *mapping = _astrometryModel->getMapping(ccdImage);
123 // count parameters
124 std::size_t npar_mapping = (_fittingDistortions) ? mapping->getNpar() : 0;
125 std::size_t npar_pos = (_fittingPos) ? 2 : 0;
126 std::size_t npar_pm = (_fittingPM) ? 2 : 0;
127 std::size_t npar_tot = npar_mapping + npar_pos + npar_pm;
128 // if (npar_tot == 0) this CcdImage does not contribute
129 // any constraint to the fit, so :
130 if (npar_tot == 0) return;
131 IndexVector indices(npar_tot, -1);
132 if (_fittingDistortions) mapping->getMappingIndices(indices);
133
134 // FittedStar is "observed" epoch, MeasuredStar is "baseline"
135 double deltaYears = _epoch - ccdImage.getEpoch();
136 // transformation from sky to TP
137 auto sky2TP = _astrometryModel->getSkyToTangentPlane(ccdImage);
138 // reserve matrices once for all measurements
139 AstrometryTransformLinear dypdy;
140 // the shape of H (et al) is required this way in order to be able to
141 // separate derivatives along x and y as vectors.
142 Eigen::MatrixX2d H(npar_tot, 2), halpha(npar_tot, 2), HW(npar_tot, 2);
143 Eigen::Matrix2d transW(2, 2);
144 Eigen::Matrix2d alpha(2, 2);
145 Eigen::VectorXd grad(npar_tot);
146 // current position in the Jacobian
147 Eigen::Index kTriplets = tripletList.getNextFreeIndex();
148 const MeasuredStarList &catalog = (msList) ? *msList : ccdImage.getCatalogForFit();
149
150 for (auto &i : catalog) {
151 const MeasuredStar &ms = *i;
152 if (!ms.isValid()) continue;
153 // tweak the measurement errors
154 FatPoint inPos = ms;
155 tweakAstromMeasurementErrors(inPos, ms, _posError);
156 H.setZero(); // we cannot be sure that all entries will be overwritten.
157 FatPoint outPos;
158 // should *not* fill H if whatToFit excludes mapping parameters.
159 if (_fittingDistortions)
160 mapping->computeTransformAndDerivatives(inPos, outPos, H);
161 else
162 mapping->transformPosAndErrors(inPos, outPos);
163
164 std::size_t ipar = npar_mapping;
165 double det = outPos.vx * outPos.vy - std::pow(outPos.vxy, 2);
166 if (det <= 0 || outPos.vx <= 0 || outPos.vy <= 0) {
167 LOGLS_WARN(_log, "Inconsistent measurement errors: drop measurement at "
168 << Point(ms) << " in image " << ccdImage.getName());
169 continue;
170 }
171 transW(0, 0) = outPos.vy / det;
172 transW(1, 1) = outPos.vx / det;
173 transW(0, 1) = transW(1, 0) = -outPos.vxy / det;
174 // compute alpha, a triangular square root
175 // of transW (i.e. a Cholesky factor)
176 alpha(0, 0) = sqrt(transW(0, 0));
177 // checked that alpha*alphaT = transW
178 alpha(1, 0) = transW(0, 1) / alpha(0, 0);
179 // DB - I think that the next line is equivalent to : alpha(1,1) = 1./sqrt(outPos.vy)
180 // PA - seems correct !
181 alpha(1, 1) = 1. / sqrt(det * transW(0, 0));
182 alpha(0, 1) = 0;
183
184 std::shared_ptr<FittedStar const> const fs = ms.getFittedStar();
185
186 Point fittedStarInTP = transformFittedStar(*fs, *sky2TP, deltaYears);
187
188 // compute derivative of TP position w.r.t sky position ....
189 if (npar_pos > 0) // ... if actually fitting FittedStar position
190 {
191 sky2TP->computeDerivative(*fs, dypdy, 1e-3);
192 // sign checked
193 // TODO Still have to check with non trivial non-diagonal terms
194 H(npar_mapping, 0) = -dypdy.A11();
195 H(npar_mapping + 1, 0) = -dypdy.A12();
196 H(npar_mapping, 1) = -dypdy.A21();
197 H(npar_mapping + 1, 1) = -dypdy.A22();
198 indices[npar_mapping] = fs->getIndexInMatrix();
199 indices.at(npar_mapping + 1) = fs->getIndexInMatrix() + 1;
200 ipar += npar_pos;
201 }
202 /* only consider proper motions of objects allowed to move,
203 unless the fit is going to be degenerate */
204 // TODO: left as reference for when we implement PM fitting
205 // TODO: mjd here would become either deltaYears or maybe associations.epoch? Check the math first!
206 // if (_fittingPM && fs->mightMove) {
207 // H(ipar, 0) = -mjd; // Sign unchecked but consistent with above
208 // H(ipar + 1, 1) = -mjd;
209 // indices[ipar] = fs->getIndexInMatrix() + 2;
210 // indices[ipar + 1] = fs->getIndexInMatrix() + 3;
211 // ipar += npar_pm;
212 // }
213
214 // We can now compute the residual
215 Eigen::Vector2d res(fittedStarInTP.x - outPos.x, fittedStarInTP.y - outPos.y);
216
217 // do not write grad = H*transW*res to avoid
218 // dynamic allocation of a temporary
219 halpha = H * alpha;
220 HW = H * transW;
221 grad = HW * res;
222 // now feed in triplets and fullGrad
223 for (std::size_t ipar = 0; ipar < npar_tot; ++ipar) {
224 for (std::size_t ic = 0; ic < 2; ++ic) {
225 double val = halpha(ipar, ic);
226 if (val == 0) continue;
227 tripletList.addTriplet(indices[ipar], kTriplets + ic, val);
228 }
229 fullGrad(indices[ipar]) += grad(ipar);
230 }
231 kTriplets += 2; // each measurement contributes 2 columns in the Jacobian
232 } // end loop on measurements
233 tripletList.setNextFreeIndex(kTriplets);
234}
T pow(T... args)
T sqrt(T... args)
ImageT val
Definition CR.cc:146

◆ leastSquareDerivativesReference()

void lsst::jointcal::AstrometryFit::leastSquareDerivativesReference ( FittedStarList const & fittedStarList,
TripletList & tripletList,
Eigen::VectorXd & grad ) const
overrideprivatevirtual

Compute the derivatives of the reference terms.

Implements lsst::jointcal::FitterBase.

Definition at line 236 of file AstrometryFit.cc.

238 {
239 /**********************************************************************/
240 /* @note the math in this method and accumulateStatRefStars() must be kept consistent,
241 * in terms of +/- convention, definition of model, etc. */
242 /**********************************************************************/
243
244 /* We compute here the derivatives of the terms involving fitted
245 stars and reference stars. They only provide contributions if we
246 are fitting positions: */
247 if (!_fittingPos) return;
248 /* the other case where the accumulation of derivatives stops
249 here is when there are no RefStars */
250 if (_associations->refStarList.size() == 0) return;
251 Eigen::Matrix2d W(2, 2);
252 Eigen::Matrix2d alpha(2, 2);
253 Eigen::Matrix2d H(2, 2), halpha(2, 2), HW(2, 2);
254 AstrometryTransformLinear der;
255 Eigen::Vector2d res, grad;
256 Eigen::Index indices[2];
257 Eigen::Index kTriplets = tripletList.getNextFreeIndex();
258 /* We cannot use the spherical coordinates directly to evaluate
259 Euclidean distances, we have to use a projector on some plane in
260 order to express least squares. Not projecting could lead to a
261 disaster around the poles or across alpha=0. So we need a
262 projector. We construct a projector and will change its
263 projection point at every object */
264 TanRaDecToPixel proj(AstrometryTransformLinear(), Point(0., 0.));
265 for (auto const &i : fittedStarList) {
266 const FittedStar &fs = *i;
267 auto rs = fs.getRefStar();
268 if (rs == nullptr) continue;
269 proj.setTangentPoint(fs);
270 // fs projects to (0,0), no need to compute its transform.
271 FatPoint rsProj;
272 proj.transformPosAndErrors(*rs, rsProj);
273 // Compute the derivative of the projector to incorporate its effects on the errors.
274 proj.computeDerivative(fs, der, 1e-4);
275 // sign checked. TODO check that the off-diagonal terms are OK.
276 H(0, 0) = -der.A11();
277 H(1, 0) = -der.A12();
278 H(0, 1) = -der.A21();
279 H(1, 1) = -der.A22();
280 // TO DO : account for proper motions.
281 double det = rsProj.vx * rsProj.vy - std::pow(rsProj.vxy, 2);
282 if (rsProj.vx <= 0 || rsProj.vy <= 0 || det <= 0) {
283 LOGLS_WARN(_log, "RefStar error matrix not positive definite for: " << *rs);
284 continue;
285 }
286 W(0, 0) = rsProj.vy / det;
287 W(0, 1) = W(1, 0) = -rsProj.vxy / det;
288 W(1, 1) = rsProj.vx / det;
289 // compute alpha, a triangular square root
290 // of W (i.e. a Cholesky factor)
291 alpha(0, 0) = sqrt(W(0, 0));
292 // checked that alpha*alphaT = transW
293 alpha(1, 0) = W(0, 1) / alpha(0, 0);
294 alpha(1, 1) = 1. / sqrt(det * W(0, 0));
295 alpha(0, 1) = 0;
296 indices[0] = fs.getIndexInMatrix();
297 indices[1] = fs.getIndexInMatrix() + 1;
298 unsigned npar_tot = 2;
299 /* TODO: account here for proper motions in the reference
300 catalog. We can code the effect and set the value to 0. Most
301 (all?) catalogs do not even come with a reference epoch. Gaia
302 will change that. When refraction enters into the game, one should
303 pay attention to the orientation of the frame */
304
305 /* The residual should be Proj(fs)-Proj(*rs) in order to be consistent
306 with the measurement terms. Since P(fs) = 0, we have: */
307 res[0] = -rsProj.x;
308 res[1] = -rsProj.y;
309 halpha = H * alpha;
310 // grad = H*W*res
311 HW = H * W;
312 grad = HW * res;
313 // now feed in triplets and fullGrad
314 for (std::size_t ipar = 0; ipar < npar_tot; ++ipar) {
315 for (unsigned ic = 0; ic < 2; ++ic) {
316 double val = halpha(ipar, ic);
317 if (val == 0) continue;
318 tripletList.addTriplet(indices[ipar], kTriplets + ic, val);
319 }
320 fullGrad(indices[ipar]) += grad(ipar);
321 }
322 kTriplets += 2; // each measurement contributes 2 columns in the Jacobian
323 }
324 tripletList.setNextFreeIndex(kTriplets);
325}

◆ minimize()

MinimizeResult lsst::jointcal::FitterBase::minimize ( std::string const & whatToFit,
double nSigmaCut = 0,
double sigmaRelativeTolerance = 0,
bool doRankUpdate = true,
bool 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]whatToFitSee child method assignIndices for valid string values.
[in]nSigmaCutHow many sigma to reject outliers at. Outlier rejection ignored for nSigmaCut=0.
[in]sigmaRelativeTolerancePercentage change in the chi2 cut for outliers tolerated for termination. If value is zero, minimization iterations will continue until there are no outliers.
[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 179 of file FitterBase.cc.

181 {
182 assignIndices(whatToFit);
183
184 MinimizeResult returnCode = MinimizeResult::Converged;
185
186 // For the initial vector size, use all measured stars + all fitted stars, which should give the
187 // maximum possible number of triplets.
188 std::size_t nTrip = (_lastNTrip)
189 ? _lastNTrip
190 : _associations->getMaxMeasuredStars() + _associations->fittedStarList.size();
191 TripletList tripletList(nTrip);
192 Eigen::VectorXd grad(_nTotal);
193 grad.setZero();
194 double scale = 1.0;
195
196 // Fill the triplets
197 leastSquareDerivatives(tripletList, grad);
198 _lastNTrip = tripletList.size();
199
200 LOGLS_DEBUG(_log, "End of triplet filling, ntrip = " << tripletList.size());
201
202 SparseMatrixD hessian = createHessian(_nTotal, tripletList);
203 tripletList.clear(); // we don't need it any more after we have the hessian.
204
205 LOGLS_DEBUG(_log, "Starting factorization, hessian: dim="
206 << hessian.rows() << " non-zeros=" << hessian.nonZeros()
207 << " filling-frac = " << hessian.nonZeros() / std::pow(hessian.rows(), 2));
208
209 if (dumpMatrixFile != "") {
210 if (hessian.rows() * hessian.cols() > 2e8) {
211 LOGLS_WARN(_log, "Hessian matrix is too big to dump to file, with rows, columns: "
212 << hessian.rows() << ", " << hessian.cols());
213 } else {
214 dumpMatrixAndGradient(hessian, grad, dumpMatrixFile, _log);
215 }
216 }
217
219 if (chol.info() != Eigen::Success) {
220 LOGLS_ERROR(_log, "minimize: factorization failed ");
222 }
223
224 std::size_t totalMeasOutliers = 0;
225 std::size_t totalRefOutliers = 0;
226 double oldChi2 = computeChi2().chi2;
227 double oldSigmaCut = 0.;
228 double sigmaCut = 0.;
229
230 while (true) {
231 Eigen::VectorXd delta = chol.solve(grad);
232 if (doLineSearch) {
233 scale = _lineSearch(delta);
234 }
235 offsetParams(scale * delta);
236 Chi2Statistic currentChi2(computeChi2());
237 LOGLS_DEBUG(_log, currentChi2);
238 if (!isfinite(currentChi2.chi2)) {
239 LOGL_ERROR(_log, "chi2 is not finite. Aborting outlier rejection.");
240 returnCode = MinimizeResult::NonFinite;
241 break;
242 }
243 if (currentChi2.chi2 > oldChi2 && totalMeasOutliers + totalRefOutliers != 0) {
244 LOGL_WARN(_log, "chi2 went up, skipping outlier rejection loop");
246 break;
247 }
248 oldChi2 = currentChi2.chi2;
249
250 if (nSigmaCut == 0) break; // no rejection step to perform
251 MeasuredStarList msOutliers;
252 FittedStarList fsOutliers;
253 // keep nOutliers so we don't have to sum msOutliers.size()+fsOutliers.size() twice below.
254 std::size_t nOutliers = findOutliers(nSigmaCut, msOutliers, fsOutliers, sigmaCut);
255 double relChange = 0.;
256 if(oldSigmaCut!=0.) relChange = (1 - sigmaCut / oldSigmaCut);
257
258 LOGLS_DEBUG(_log, "findOutliers chi2 cut level: " << sigmaCut << ", relative change: " << relChange);
259 // If sigmaRelativeTolerance is set and at least one iteration has been done, break loop when the
260 // fractional change in sigmaCut levels is less than the sigmaRelativeTolerance parameter.
261 if ((sigmaRelativeTolerance > 0) && (oldSigmaCut > 0 && relChange < sigmaRelativeTolerance)) {
262 LOGLS_INFO(_log, "Iterations stopped with chi2 cut at " << sigmaCut << " and relative change of "
263 << relChange);
264 break;
265 }
266 totalMeasOutliers += msOutliers.size();
267 totalRefOutliers += fsOutliers.size();
268 oldSigmaCut = sigmaCut;
269 if (nOutliers == 0) break;
270 TripletList outlierTriplets(nOutliers);
271 grad.setZero(); // recycle the gradient
272 // compute the contributions of outliers to derivatives
273 outliersContributions(msOutliers, fsOutliers, outlierTriplets, grad);
274 // Remove significant outliers
275 removeMeasOutliers(msOutliers);
276 removeRefOutliers(fsOutliers);
277 if (doRankUpdate) {
278 // convert triplet list to eigen internal format
279 SparseMatrixD H(_nTotal, outlierTriplets.getNextFreeIndex());
280 H.setFromTriplets(outlierTriplets.begin(), outlierTriplets.end());
281 chol.update(H, false /* means downdate */);
282 // The contribution of outliers to the gradient is the opposite
283 // of the contribution of all other terms, because they add up to 0
284 grad *= -1;
285 } else {
286 // don't reuse tripletList because we want a new nextFreeIndex.
287 TripletList nextTripletList(_lastNTrip);
288 grad.setZero();
289 // Rebuild the matrix and gradient
290 leastSquareDerivatives(nextTripletList, grad);
291 _lastNTrip = nextTripletList.size();
292 LOGLS_DEBUG(_log, "Triplets recomputed, ntrip = " << nextTripletList.size());
293
294 hessian = createHessian(_nTotal, nextTripletList);
295 nextTripletList.clear(); // we don't need it any more after we have the hessian.
296
298 "Restarting factorization, hessian: dim="
299 << hessian.rows() << " non-zeros=" << hessian.nonZeros()
300 << " filling-frac = " << hessian.nonZeros() / std::pow(hessian.rows(), 2));
301 chol.compute(hessian);
302 if (chol.info() != Eigen::Success) {
303 LOGLS_ERROR(_log, "minimize: factorization failed ");
305 }
306 }
307 }
308
309 if (totalMeasOutliers + totalRefOutliers > 0) {
310 _associations->cleanFittedStars();
311 }
312
313 // only print the outlier summary if outlier rejection was turned on.
314 if (nSigmaCut != 0) {
315 LOGLS_INFO(_log, "Number of outliers (Measured + Reference = Total): "
316 << totalMeasOutliers << " + " << totalRefOutliers << " = "
317 << totalMeasOutliers + totalRefOutliers);
318 }
319 return returnCode;
320}
#define LOGL_WARN(logger, message...)
Log a warn-level message using a varargs/printf style interface.
Definition Log.h:547
#define LOGLS_ERROR(logger, message)
Log a error-level message using an iostream-based interface.
Definition Log.h:679
#define LOGL_ERROR(logger, message...)
Log a error-level message using a varargs/printf style interface.
Definition Log.h:563
void removeRefOutliers(FittedStarList &outliers)
Remove refStar outliers from the fit. No Refit done.
Chi2Statistic computeChi2() const
Returns the chi2 for the current state.
Definition FitterBase.cc:43
virtual void assignIndices(std::string const &whatToFit)=0
Set parameters to fit and assign indices in the big matrix.
virtual void offsetParams(Eigen::VectorXd const &delta)=0
Offset the parameters by the requested quantities.
void outliersContributions(MeasuredStarList &msOutliers, FittedStarList &fsOutliers, TripletList &tripletList, Eigen::VectorXd &grad)
Contributions to derivatives from (presumably) outlier 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.
Definition FitterBase.cc:53
void removeMeasOutliers(MeasuredStarList &outliers)
Remove measuredStar outliers from the fit. No Refit done.
T isfinite(T... args)
scale(algorithm, min, max=None, frame=None)
Definition ds9.py:108

◆ offsetParams()

void lsst::jointcal::AstrometryFit::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]deltavector of offsets to apply

Implements lsst::jointcal::FitterBase.

Definition at line 454 of file AstrometryFit.cc.

454 {
455 if (delta.size() != _nTotal)
456 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
457 "AstrometryFit::offsetParams : the provided vector length is not compatible with "
458 "the current whatToFit setting");
459 if (_fittingDistortions) _astrometryModel->offsetParams(delta);
460
461 if (_fittingPos) {
462 FittedStarList &fittedStarList = _associations->fittedStarList;
463 for (auto const &i : fittedStarList) {
464 FittedStar &fs = *i;
465 // the parameter layout here is used also
466 // - when filling the derivatives
467 // - when assigning indices (assignIndices())
468 Eigen::Index index = fs.getIndexInMatrix();
469 fs.x += delta(index);
470 fs.y += delta(index + 1);
471 // TODO: left as reference for when we implement PM fitting
472 // if ((_fittingPM)&fs.mightMove) {
473 // fs.pmx += delta(index + 2);
474 // fs.pmy += delta(index + 3);
475 // }
476 }
477 }
478}
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition Exception.h:48

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ outliersContributions()

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

Contributions to derivatives from (presumably) outlier terms.

No discarding done.

Definition at line 322 of file FitterBase.cc.

323 {
324 for (auto &outlier : msOutliers) {
325 MeasuredStarList tmp;
326 tmp.push_back(outlier);
327 const CcdImage &ccdImage = outlier->getCcdImage();
328 leastSquareDerivativesMeasurement(ccdImage, tripletList, grad, &tmp);
329 }
330 leastSquareDerivativesReference(fsOutliers, tripletList, grad);
331}

◆ removeMeasOutliers()

void lsst::jointcal::FitterBase::removeMeasOutliers ( MeasuredStarList & outliers)
protectedinherited

Remove measuredStar outliers from the fit. No Refit done.

Definition at line 333 of file FitterBase.cc.

333 {
334 for (auto &measuredStar : outliers) {
335 auto fittedStar = measuredStar->getFittedStar();
336 measuredStar->setValid(false);
337 fittedStar->getMeasurementCount()--; // could be put in setValid
338 }
339}

◆ removeRefOutliers()

void lsst::jointcal::FitterBase::removeRefOutliers ( FittedStarList & outliers)
protectedinherited

Remove refStar outliers from the fit. No Refit done.

Definition at line 341 of file FitterBase.cc.

341 {
342 for (auto &fittedStar : outliers) {
343 fittedStar->setRefStar(nullptr);
344 }
345}

◆ 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 355 of file FitterBase.cc.

355 {
356 std::string replaceStr = "{type}";
357 auto pos = baseName.find(replaceStr);
358 std::string measFilename(baseName);
359 measFilename.replace(pos, replaceStr.size(), "-meas.csv");
360 std::string refFilename(baseName);
361 refFilename.replace(pos, replaceStr.size(), "-ref.csv");
362 saveChi2MeasContributions(measFilename);
363 saveChi2RefContributions(refFilename);
364}
virtual void saveChi2MeasContributions(std::string const &filename) const =0
Save a CSV file containing residuals of measurement terms.
virtual void saveChi2RefContributions(std::string const &filename) const =0
Save a CSV file containing residuals of reference terms.
T size(T... args)

◆ saveChi2MeasContributions()

void lsst::jointcal::AstrometryFit::saveChi2MeasContributions ( std::string const & filename) const
overrideprotectedvirtual

Save a CSV file containing residuals of measurement terms.

Implements lsst::jointcal::FitterBase.

Definition at line 496 of file AstrometryFit.cc.

496 {
497 std::ofstream ofile(filename.c_str());
498 std::string separator = "\t";
499
500 ofile << "id" << separator << "xccd" << separator << "yccd " << separator;
501 ofile << "rx" << separator << "ry" << separator;
502 ofile << "xtp" << separator << "ytp" << separator;
503 ofile << "mag" << separator << "deltaYears" << separator;
504 ofile << "xErr" << separator << "yErr" << separator << "xyCov" << separator;
505 ofile << "xtpi" << separator << "ytpi" << separator;
506 ofile << "rxi" << separator << "ryi" << separator;
507 ofile << "fsindex" << separator;
508 ofile << "ra" << separator << "dec" << separator;
509 ofile << "chi2" << separator << "nm" << separator;
510 ofile << "detector" << separator << "visit" << std::endl;
511
512 ofile << "id in source catalog" << separator << "coordinates in CCD (pixels)" << separator << separator;
513 ofile << "residual on TP (degrees)" << separator << separator;
514 ofile << "transformed coordinate in TP (degrees)" << separator << separator;
515 ofile << "rough magnitude" << separator << "Julian epoch year delta from fit epoch" << separator;
516 ofile << "transformed measurement uncertainty (degrees)" << separator << separator << separator;
517 ofile << "as-read position on TP (degrees)" << separator << separator;
518 ofile << "as-read residual on TP (degrees)" << separator << separator;
519 ofile << "unique index of the fittedStar" << separator;
520 ofile << "on sky position of fittedStar" << separator << separator;
521 ofile << "contribution to chi2 (2D dofs)" << separator << "number of measurements of this fittedStar"
522 << separator;
523 ofile << "detector id" << separator << "visit id" << std::endl;
524
525 const CcdImageList &ccdImageList = _associations->getCcdImageList();
526 for (auto const &ccdImage : ccdImageList) {
527 const MeasuredStarList &cat = ccdImage->getCatalogForFit();
528 const AstrometryMapping *mapping = _astrometryModel->getMapping(*ccdImage);
529 const auto readTransform = ccdImage->getReadWcs();
530 // FittedStar is "observed" epoch, MeasuredStar is "baseline"
531 double deltaYears = _epoch - ccdImage->getEpoch();
532 for (auto const &ms : cat) {
533 if (!ms->isValid()) continue;
534 FatPoint tpPos;
535 FatPoint inPos = *ms;
536 tweakAstromMeasurementErrors(inPos, *ms, _posError);
537 mapping->transformPosAndErrors(inPos, tpPos);
538 auto sky2TP = _astrometryModel->getSkyToTangentPlane(*ccdImage);
539 const std::unique_ptr<AstrometryTransform> readPixToTangentPlane =
540 compose(*sky2TP, *readTransform);
541 FatPoint inputTpPos = readPixToTangentPlane->apply(inPos);
542 std::shared_ptr<FittedStar const> const fs = ms->getFittedStar();
543
544 Point fittedStarInTP = transformFittedStar(*fs, *sky2TP, deltaYears);
545 Point res = tpPos - fittedStarInTP;
546 Point inputRes = inputTpPos - fittedStarInTP;
547 double det = tpPos.vx * tpPos.vy - std::pow(tpPos.vxy, 2);
548 double wxx = tpPos.vy / det;
549 double wyy = tpPos.vx / det;
550 double wxy = -tpPos.vxy / det;
551 double chi2 = wxx * res.x * res.x + wyy * res.y * res.y + 2 * wxy * res.x * res.y;
552
553 ofile << std::setprecision(9);
554 ofile << ms->getId() << separator << ms->x << separator << ms->y << separator;
555 ofile << res.x << separator << res.y << separator;
556 ofile << tpPos.x << separator << tpPos.y << separator;
557 ofile << fs->getMag() << separator << deltaYears << separator;
558 ofile << tpPos.vx << separator << tpPos.vy << separator << tpPos.vxy << separator;
559 ofile << inputTpPos.x << separator << inputTpPos.y << separator;
560 ofile << inputRes.x << separator << inputRes.y << separator;
561 ofile << fs->getIndexInMatrix() << separator;
562 ofile << fs->x << separator << fs->y << separator;
563 ofile << chi2 << separator << fs->getMeasurementCount() << separator;
564 ofile << ccdImage->getCcdId() << separator << ccdImage->getVisit() << std::endl;
565 } // loop on measurements in image
566 } // loop on images
567}
T endl(T... args)
std::unique_ptr< AstrometryTransform > compose(AstrometryTransform const &left, AstrometryTransform const &right)
Returns a pointer to a composition of transforms, representing left(right()).
std::list< std::shared_ptr< CcdImage > > CcdImageList
Definition CcdImage.h:46
T setprecision(T... args)

◆ saveChi2RefContributions()

void lsst::jointcal::AstrometryFit::saveChi2RefContributions ( std::string const & filename) const
overrideprotectedvirtual

Save a CSV file containing residuals of reference terms.

Implements lsst::jointcal::FitterBase.

Definition at line 569 of file AstrometryFit.cc.

569 {
570 std::ofstream ofile(filename.c_str());
571 std::string separator = "\t";
572
573 ofile << "ra" << separator << "dec " << separator;
574 ofile << "rx" << separator << "ry" << separator;
575 ofile << "mag" << separator;
576 ofile << "xErr" << separator << "yErr" << separator << "xyCov" << separator;
577 ofile << "fsindex" << separator;
578 ofile << "chi2" << separator << "nm" << std::endl;
579
580 ofile << "coordinates of fittedStar (degrees)" << separator << separator;
581 ofile << "residual on TP (degrees)" << separator << separator;
582 ofile << "magnitude" << separator;
583 ofile << "refStar transformed measurement uncertainty (degrees)" << separator << separator << separator;
584 ofile << "unique index of the fittedStar" << separator;
585 ofile << "refStar contribution to chi2 (2D dofs)" << separator
586 << "number of measurements of this FittedStar" << std::endl;
587
588 // The following loop is heavily inspired from AstrometryFit::computeChi2()
589 const FittedStarList &fittedStarList = _associations->fittedStarList;
590 TanRaDecToPixel proj(AstrometryTransformLinear(), Point(0., 0.));
591 for (auto const &i : fittedStarList) {
592 const FittedStar &fs = *i;
593 auto rs = fs.getRefStar();
594 if (rs == nullptr) continue;
595 proj.setTangentPoint(fs);
596 // fs projects to (0,0), no need to compute its transform.
597 FatPoint rsProj;
598 proj.transformPosAndErrors(*rs, rsProj);
599 double chi2 = computeProjectedRefStarChi2(rsProj);
600
601 ofile << std::setprecision(9);
602 ofile << fs.x << separator << fs.y << separator;
603 ofile << rsProj.x << separator << rsProj.y << separator;
604 ofile << fs.getMag() << separator;
605 ofile << rsProj.vx << separator << rsProj.vy << separator << rsProj.vxy << separator;
606 ofile << fs.getIndexInMatrix() << separator;
607 ofile << chi2 << separator << fs.getMeasurementCount() << std::endl;
608 } // loop on FittedStars
609}

Member Data Documentation

◆ _associations

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

Definition at line 165 of file FitterBase.h.

◆ _lastNTrip

Eigen::Index lsst::jointcal::FitterBase::_lastNTrip
protectedinherited

Definition at line 168 of file FitterBase.h.

◆ _log

LOG_LOGGER lsst::jointcal::FitterBase::_log
protectedinherited

Definition at line 174 of file FitterBase.h.

◆ _nModelParams

Eigen::Index lsst::jointcal::FitterBase::_nModelParams
protectedinherited

Definition at line 170 of file FitterBase.h.

◆ _nStarParams

Eigen::Index lsst::jointcal::FitterBase::_nStarParams
protectedinherited

Definition at line 171 of file FitterBase.h.

◆ _nTotal

Eigen::Index lsst::jointcal::FitterBase::_nTotal
protectedinherited

Definition at line 169 of file FitterBase.h.

◆ _whatToFit

std::string lsst::jointcal::FitterBase::_whatToFit
protectedinherited

Definition at line 166 of file FitterBase.h.


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