lsst.jointcal  14.0-17-ge3cc87b
PhotometryFit.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <iomanip>
3 #include <algorithm>
4 #include <fstream>
5 
6 #include "Eigen/Sparse"
7 
8 #include "lsst/log/Log.h"
9 #include "lsst/pex/exceptions.h"
12 #include "lsst/jointcal/Chi2.h"
14 #include "lsst/jointcal/Gtransfo.h"
16 
17 namespace {
18 LOG_LOGGER _log = LOG_GET("jointcal.PhotometryFit");
19 }
20 
21 namespace lsst {
22 namespace jointcal {
23 
24 void PhotometryFit::leastSquareDerivativesMeasurement(CcdImage const &ccdImage, TripletList &tripletList,
25  Eigen::VectorXd &grad,
26  MeasuredStarList const *measuredStarList) const {
27  /**********************************************************************/
28  /* @note the math in this method and accumulateStatImageList() must be kept consistent,
29  * in terms of +/- convention, definition of model, etc. */
30  /**********************************************************************/
31 
32  /* this routine works in two different ways: either providing the
33  Ccd, of providing the MeasuredStarList. In the latter case, the
34  Ccd should match the one(s) in the list. */
35  if (measuredStarList) assert(&(measuredStarList->front()->getCcdImage()) == &ccdImage);
36 
37  unsigned nparModel = (_fittingModel) ? _photometryModel->getNpar(ccdImage) : 0;
38  unsigned nparFlux = (_fittingFluxes) ? 1 : 0;
39  unsigned nparTotal = nparModel + nparFlux;
40  std::vector<unsigned> indices(nparModel, -1);
41  if (_fittingModel) _photometryModel->getMappingIndices(ccdImage, indices);
42 
43  Eigen::VectorXd H(nparTotal); // derivative matrix
44  // current position in the Jacobian
45  unsigned kTriplets = tripletList.getNextFreeIndex();
46  const MeasuredStarList &catalog = (measuredStarList) ? *measuredStarList : ccdImage.getCatalogForFit();
47 
48  for (auto const &measuredStar : catalog) {
49  if (!measuredStar->isValid()) continue;
50 // tweak the measurement errors
51 #ifdef FUTURE
52  TweakPhotomMeasurementErrors(inPos, *measuredStar, _fluxError);
53 #endif
54  H.setZero(); // we cannot be sure that all entries will be overwritten.
55 
56  double residual = _photometryModel->transform(ccdImage, *measuredStar, measuredStar->getInstFlux()) -
57  measuredStar->getFittedStar()->getFlux();
58 
59  double inverseSigma =
60  1.0 / _photometryModel->transform(ccdImage, *measuredStar, measuredStar->getInstFluxErr());
61  double W = std::pow(inverseSigma, 2);
62 
63  if (_fittingModel) {
64  _photometryModel->computeParameterDerivatives(*measuredStar, ccdImage, H);
65  for (unsigned k = 0; k < indices.size(); k++) {
66  unsigned l = indices[k];
67  tripletList.addTriplet(l, kTriplets, H[k] * inverseSigma);
68  grad[l] += H[k] * W * residual;
69  }
70  }
71  if (_fittingFluxes) {
72  unsigned index = measuredStar->getFittedStar()->getIndexInMatrix();
73  // Note: H = dR/dFittedStarFlux == -1
74  tripletList.addTriplet(index, kTriplets, -1.0 * inverseSigma);
75  grad[index] += -1.0 * W * residual;
76  }
77  kTriplets += 1; // each measurement contributes 1 column in the Jacobian
78  }
79 
80  tripletList.setNextFreeIndex(kTriplets);
81 }
82 
83 void PhotometryFit::leastSquareDerivativesReference(FittedStarList const &fittedStarList,
84  TripletList &tripletList, Eigen::VectorXd &grad) const {
85  /**********************************************************************/
88  /**********************************************************************/
89 
90  // Derivatives of terms involving fitted and refstars only contribute if we are fitting fluxes.
91  if (!_fittingFluxes) return;
92  // Can't compute anything if there are no refStars.
93  if (_associations->refStarList.size() == 0) return;
94 
95  unsigned kTriplets = tripletList.getNextFreeIndex();
96 
97  for (auto const &fittedStar : fittedStarList) {
98  auto refStar = fittedStar->getRefStar();
99  if (refStar == nullptr) continue; // no contribution if no associated refstar
100  // TODO: Can we actually work with multiple filters at a time in this scheme?
101  // TODO: I feel we might only be able to fit one filter at a time, because these terms are
102  // independent of the ccdImage, so we don't have a specific filter.
103  // filter = ccdImage.getFilter();
104 
105  // W == inverseSigma^2
106  double inverseSigma = 1.0 / refStar->getFluxErr();
107  // Residual is fittedStar.flux - refStar.flux for consistency with measurement terms.
108  double residual = fittedStar->getFlux() - refStar->getFlux();
109 
110  unsigned index = fittedStar->getIndexInMatrix();
111  // Note: H = dR/dFittedStar == 1
112  tripletList.addTriplet(index, kTriplets, 1.0 * inverseSigma);
113  grad(index) += 1.0 * std::pow(inverseSigma, 2) * residual;
114  kTriplets += 1;
115  }
116  tripletList.setNextFreeIndex(kTriplets);
117 }
118 
119 void PhotometryFit::accumulateStatImageList(CcdImageList const &ccdImageList, Chi2Accumulator &accum) const {
120  /**********************************************************************/
123  /**********************************************************************/
124  for (auto const &ccdImage : ccdImageList) {
125  auto &catalog = ccdImage->getCatalogForFit();
126 
127  for (auto const &measuredStar : catalog) {
128  if (!measuredStar->isValid()) continue;
129  double sigma =
130  _photometryModel->transform(*ccdImage, *measuredStar, measuredStar->getInstFluxErr());
131 #ifdef FUTURE
132  TweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
133 #endif
134  double residual =
135  _photometryModel->transform(*ccdImage, *measuredStar, measuredStar->getInstFlux()) -
136  measuredStar->getFittedStar()->getFlux();
137 
138  double chi2Val = std::pow(residual / sigma, 2);
139  accum.addEntry(chi2Val, 1, measuredStar);
140  } // end loop on measurements
141  }
142 }
143 
144 void PhotometryFit::accumulateStatRefStars(Chi2Accumulator &accum) const {
145  /**********************************************************************/
148  /**********************************************************************/
149 
150  FittedStarList &fittedStarList = _associations->fittedStarList;
151  for (auto const &fittedStar : fittedStarList) {
152  auto refStar = fittedStar->getRefStar();
153  if (refStar == nullptr) continue;
154  double chi2 = std::pow(((fittedStar->getFlux() - refStar->getFlux()) / refStar->getFluxErr()), 2);
155  accum.addEntry(chi2, 1, fittedStar);
156  }
157 }
158 
160 
162 void PhotometryFit::getIndicesOfMeasuredStar(MeasuredStar const &measuredStar,
163  std::vector<unsigned> &indices) const {
164  indices.clear();
165  if (_fittingModel) {
166  _photometryModel->getMappingIndices(measuredStar.getCcdImage(), indices);
167  }
168  if (_fittingFluxes) {
169  auto fs = measuredStar.getFittedStar();
170  unsigned fsIndex = fs->getIndexInMatrix();
171  indices.push_back(fsIndex);
172  }
173 }
174 
176  _whatToFit = whatToFit;
177  LOGLS_INFO(_log, "assignIndices: now fitting: " << whatToFit);
178  _fittingModel = (_whatToFit.find("Model") != std::string::npos);
179  _fittingFluxes = (_whatToFit.find("Fluxes") != std::string::npos);
180  // When entering here, we assume that whatToFit has already been interpreted.
181 
182  _nParModel = (_fittingModel) ? _photometryModel->assignIndices(whatToFit, 0) : 0;
183  unsigned ipar = _nParModel;
184 
185  if (_fittingFluxes) {
186  for (auto &fittedStar : _associations->fittedStarList) {
187  // the parameter layout here is used also
188  // - when filling the derivatives
189  // - when updating (offsetParams())
190  // - in getIndicesOfMeasuredStar
191  fittedStar->setIndexInMatrix(ipar);
192  ipar += 1;
193  }
194  }
195  _nParFluxes = ipar - _nParModel;
196  _nParTot = ipar;
197  LOGLS_DEBUG(_log,
198  "nParameters total: " << _nParTot << " model: " << _nParModel << " fluxes: " << _nParFluxes);
199 }
200 
201 void PhotometryFit::offsetParams(Eigen::VectorXd const &delta) {
202  if (delta.size() != _nParTot)
204  "PhotometryFit::offsetParams : the provided vector length is not compatible with "
205  "the current whatToFit setting");
206  if (_fittingModel) _photometryModel->offsetParams(delta);
207 
208  if (_fittingFluxes) {
209  for (auto &fittedStar : _associations->fittedStarList) {
210  // the parameter layout here is used also
211  // - when filling the derivatives
212  // - when assigning indices (assignIndices())
213  unsigned index = fittedStar->getIndexInMatrix();
214  fittedStar->getFlux() -= delta(index);
215  }
216  }
217 }
218 
220  std::ofstream ofile(baseName.c_str());
221  /* If we think the some coordinate on the focal plane is relevant in
222  the ntuple, because thmodel relies on it, then we have to add
223  some function to the model that returns this relevant
224  coordinate. */
225  std::string separator = "\t";
226  ofile << "#xccd" << separator << "yccd" << separator << "mag" << separator << "instFlux" << separator
227  << "instFluxError" << separator << "flux" << separator << "fluxError" << separator
228  << "transformedFlux" << separator << "transformedFluxErr" << separator << "fflux" << separator
229  << "mjd" << separator << "color" << separator << "fsindex" << separator << "ra" << separator
230  << "dec" << separator << "chi2" << separator << "nm" << separator << "chip" << separator << "visit"
231  << separator << std::endl;
232  ofile << "#coordinates in CCD" << separator << separator << "rough mag" << separator
233  << "measured instrument flux" << separator << "measured instrument flux error" << separator
234  << "measured flux" << separator << "measured flux error" << separator << separator << separator
235  << "fitted flux" << separator << "modified Julian date of the measurement" << separator
236  << "currently unused" << separator << "unique index of the fittedStar" << separator
237  << "pos of fitted star" << separator << separator << "contribution to Chi2 (1 dof)" << separator
238  << "number of measurements of this FittedStar" << separator << "chip id" << separator << "visit id"
239  << std::endl;
240  const CcdImageList &ccdImageList = _associations->getCcdImageList();
241  for (auto const &ccdImage : ccdImageList) {
242  const MeasuredStarList &cat = ccdImage->getCatalogForFit();
243  for (auto const &measuredStar : cat) {
244  if (!measuredStar->isValid()) continue;
245  double sigma =
246  _photometryModel->transform(*ccdImage, *measuredStar, measuredStar->getInstFluxErr());
247 #ifdef FUTURE
248  tweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
249 #endif
250  double flux = _photometryModel->transform(*ccdImage, *measuredStar, measuredStar->getInstFlux());
251  double fluxErr =
252  _photometryModel->transform(*ccdImage, *measuredStar, measuredStar->getInstFluxErr());
253  double jd = ccdImage->getMjd();
254  auto fittedStar = measuredStar->getFittedStar();
255  double residual = flux - fittedStar->getFlux();
256  double chi2Val = std::pow(residual / sigma, 2);
257 
258  ofile << measuredStar->x << separator << measuredStar->y << separator << fittedStar->getMag()
259  << separator << measuredStar->getInstFlux() << separator << measuredStar->getInstFluxErr()
260  << separator << measuredStar->getFlux() << separator << measuredStar->getFluxErr()
261  << separator << flux << separator << fluxErr << separator << fittedStar->getFlux()
262  << separator << jd << separator << fittedStar->color << separator
263  << fittedStar->getIndexInMatrix() << separator << fittedStar->x << separator
264  << fittedStar->y << separator << chi2Val << separator << fittedStar->getMeasurementCount()
265  << separator << ccdImage->getCcdId() << separator << ccdImage->getVisit() << std::endl;
266  } // loop on measurements in image
267  } // loop on images
268 }
269 
271  std::ofstream ofile(baseName.c_str());
272  std::string separator = "\t";
273  ofile << "#ra" << separator << "dec " << separator << "mag" << separator << "color" << separator
274  << "refFlux" << separator << "refFluxErr" << separator << "fittedFlux" << separator
275  << "fittedFluxErr" << separator << "fsindex" << separator << "chi2" << separator << "nm"
276  << std::endl;
277  ofile << "#coordinates of fittedStar" << separator << separator << "magnitude" << separator
278  << "currently unused" << separator << "default refStar flux" << separator
279  << "default refStar fluxErr" << separator << "fittedStar flux" << separator << "fittedStar fluxErr"
280  << separator << "unique index of the fittedStar" << separator
281  << "refStar contribution to Chi2 (2D dofs)" << separator
282  << "number of measurements of this FittedStar" << std::endl;
283  // The following loop is heavily inspired from AstrometryFit::computeChi2()
284  const FittedStarList &fittedStarList = _associations->fittedStarList;
285  for (auto const &fittedStar : fittedStarList) {
286  const RefStar *refStar = fittedStar->getRefStar();
287  if (refStar == nullptr) continue;
288 
289  double chi2 = std::pow(((fittedStar->getFlux() - refStar->getFlux()) / refStar->getFluxErr()), 2);
290 
291  ofile << std::setprecision(9);
292  ofile << fittedStar->x << separator << fittedStar->y << separator << fittedStar->getMag() << separator
293  << fittedStar->color << separator << refStar->getFlux() << separator << refStar->getFluxErr()
294  << separator << fittedStar->getFlux() << separator << fittedStar->getFluxErr() << separator
295  << fittedStar->getIndexInMatrix() << separator << chi2 << separator
296  << fittedStar->getMeasurementCount() << std::endl;
297  } // loop on FittedStars
298 }
299 
300 } // namespace jointcal
301 } // namespace lsst
Objects used as position anchors, typically USNO stars.
Definition: RefStar.h:16
double getFluxErr(size_t filter) const
reference fluxErr in a given filter
Definition: RefStar.h:44
void saveChi2RefContributions(std::string const &baseName) const override
Save a CSV file containing residuals of reference terms.
T endl(T... args)
A list of MeasuredStar. They are usually filled in Associations::AddImage.
Definition: MeasuredStar.h:111
afw::table::Key< double > sigma
#define LOGLS_INFO(logger, message)
STL class.
T push_back(T... args)
STL class.
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the requested quantities.
Class for a simple mapping implementing a generic Gtransfo.
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition: FittedStar.h:121
T clear(T... args)
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:113
void assignIndices(std::string const &whatToFit) override
Set parameters to fit and assign indices in the big matrix.
T find(T... args)
#define LSST_EXCEPT(type,...)
T pow(T... args)
T c_str(T... args)
#define LOGLS_DEBUG(logger, message)
double getFlux(size_t filter) const
reference flux in a given filter
Definition: RefStar.h:42
std::list< std::shared_ptr< CcdImage > > CcdImageList
Definition: CcdImage.h:23
T setprecision(T... args)
#define LOG_GET(logger)
void saveChi2MeasContributions(std::string const &baseName) const override
Save a CSV file containing residuals of measurement terms.