Class FitterBase

Inheritance Relationships

Derived Types

Class Documentation

class FitterBase

Base class for fitters.

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

Subclassed by lsst::jointcal::AstrometryFit, lsst::jointcal::PhotometryFit

Public Functions

FitterBase(std::shared_ptr<Associations> associations)
FitterBase(FitterBase const&)

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

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

Does a 1 step minimization, assuming a linear model.

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.

Return

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.

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)
    

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.

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.

virtual void offsetParams(Eigen::VectorXd const &delta) = 0

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

virtual void assignIndices(std::string const &whatToFit) = 0

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

Parameters
  • [in] whatToFit: See child class documentation for valid string values.

virtual void saveChi2Contributions(std::string const &baseName) const

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).

Protected Functions

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.

std::size_t findOutliers(double nSigmaCut, MeasuredStarList &msOutliers, FittedStarList &fsOutliers) const

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

\[ <chi2> + nSigmaCut + rms(chi2). \]
The outliers are NOT removed, and no refit is done.

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.

Return

Total number of outliers that were removed.

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

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

Contributions to derivatives from (presumably) outlier terms. No discarding done.

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.

virtual void getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, IndexVector &indices) const = 0

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

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.

virtual void accumulateStatRefStars(Chi2Accumulator &accum) const = 0

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

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.

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

virtual void leastSquareDerivativesReference(FittedStarList const &fittedStarList, TripletList &tripletList, Eigen::VectorXd &grad) const = 0

Compute the derivatives of the reference terms.

Protected Attributes

std::shared_ptr<Associations> _associations
std::string _whatToFit
Eigen::Index _lastNTrip
Eigen::Index _nParTot
Eigen::Index _nMeasuredStars
LOG_LOGGER _log