lsst.jointcal  14.0-12-gc7bdbdc
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 
219 void PhotometryFit::saveResultTuples(std::string const &tupleName) const {
220  std::ofstream ofile(tupleName.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  ofile << "#xccd: coordinate in CCD"
226  << "\t"
227  << "yccd: "
228  << "\t"
229  << "mag: rough mag"
230  << "\t"
231  << "instFlux : measured instrument flux"
232  << "\t"
233  << "instFluxError : measured instrument flux error"
234  << "\t"
235  << "flux : measured flux"
236  << "\t"
237  << "fluxError : measured flux error"
238  << "\t"
239  << "transformedFlux:"
240  << "\t"
241  << "transformedFluxErr:"
242  << "\t"
243  << "fflux : fitted flux"
244  << "\t"
245  << "mjd: Julian date of the measurement"
246  << "\t"
247  << "color : "
248  << "\t"
249  << "fsindex: some unique index of the object"
250  << "\t"
251  << "ra: pos of fitted star"
252  << "\t"
253  << "dec: pos of fitted star"
254  << "\t"
255  << "chi2: contribution to Chi2 (1 dof)"
256  << "\t"
257  << "nm: number of measurements of this FittedStar"
258  << "\t"
259  << "chip: chip number"
260  << "\t"
261  << "visit: visit id"
262  << "\t"
263  << "end" << std::endl;
264  const CcdImageList &ccdImageList = _associations->getCcdImageList();
265  for (auto const &i : ccdImageList) {
266  const CcdImage &ccdImage = *i;
267  const MeasuredStarList &cat = ccdImage.getCatalogForFit();
268  for (auto const &is : cat) {
269  const MeasuredStar &measuredStar = *is;
270  if (!measuredStar.isValid()) continue;
271  double sigma = measuredStar.getFluxErr();
272 #ifdef FUTURE
273  tweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
274 #endif
275  double flux = _photometryModel->transform(ccdImage, measuredStar, measuredStar.getInstFlux());
276  double fluxErr =
277  _photometryModel->transform(ccdImage, measuredStar, measuredStar.getInstFluxErr());
278  double jd = ccdImage.getMjd();
279  auto fs = measuredStar.getFittedStar();
280  double residual = flux - fs->getFlux();
281  double chi2Val = std::pow(residual / sigma, 2);
282  ofile << measuredStar.x << "\t" << measuredStar.y << "\t" << fs->getMag() << "\t"
283  << measuredStar.getInstFlux() << "\t" << measuredStar.getInstFluxErr() << "\t"
284  << measuredStar.getFlux() << "\t" << measuredStar.getFluxErr() << "\t" << flux << "\t"
285  << fluxErr << "\t" << fs->getFlux() << "\t" << jd << "\t" << fs->color << "\t"
286  << fs->getIndexInMatrix() << "\t" << fs->x << "\t" << fs->y << "\t" << chi2Val << "\t"
287  << fs->getMeasurementCount() << "\t" << ccdImage.getCcdId() << "\t" << ccdImage.getVisit()
288  << std::endl;
289  } // loop on measurements in image
290  } // loop on images
291 }
292 } // namespace jointcal
293 } // namespace lsst
VisitIdType getVisit() const
returns visit ID
Definition: CcdImage.h:101
double getMjd() const
Julian Date.
Definition: CcdImage.h:111
T endl(T... args)
CcdIdType getCcdId() const
returns ccd ID
Definition: CcdImage.h:98
double getFluxErr() const
Definition: BaseStar.h:64
A list of MeasuredStar. They are usually filled in Associations::AddImage.
Definition: MeasuredStar.h:111
afw::table::Key< double > sigma
STL class.
double x
coordinate
Definition: Point.h:18
T push_back(T... args)
STL class.
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the requested quantities.
#define LOGLS_DEBUG(logger, message)
Class for a simple mapping implementing a generic Gtransfo.
void saveResultTuples(std::string const &tupleName) const override
Save the full chi2 term per star that was used in the minimization, for debugging.
double getInstFluxErr() const
Definition: MeasuredStar.h:69
objects measured on actual images.
Definition: MeasuredStar.h:18
T clear(T... args)
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:107
double getFlux() const
Definition: BaseStar.h:60
#define LOGLS_INFO(logger, message)
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,...)
std::shared_ptr< const FittedStar > getFittedStar() const
Definition: MeasuredStar.h:83
T pow(T... args)
MeasuredStarList const & getCatalogForFit() const
Gets the catalog to be used for fitting, which may have been cleaned-up.
Definition: CcdImage.h:64
T c_str(T... args)
std::list< std::shared_ptr< CcdImage > > CcdImageList
Definition: CcdImage.h:23
Handler of an actual image from a single CCD.
Definition: CcdImage.h:34
#define LOG_GET(logger)
bool isValid() const
Fits may use that to discard outliers.
Definition: MeasuredStar.h:90