lsst.jointcal  master-ga8493ae4fe+6
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 static double sqr(double x) { return x * x; }
18 
19 namespace {
20 LOG_LOGGER _log = LOG_GET("jointcal.PhotometryFit");
21 }
22 
23 namespace {
25 double computeMeasurementResidual(double photomFactor, lsst::jointcal::MeasuredStar const &measuredStar) {
26  return photomFactor * measuredStar.getFlux() - measuredStar.getFittedStar()->getFlux();
27 }
28 } // namespace
29 
30 namespace lsst {
31 namespace jointcal {
32 
33 void PhotometryFit::leastSquareDerivativesMeasurement(CcdImage const &ccdImage, TripletList &tripletList,
34  Eigen::VectorXd &grad,
35  MeasuredStarList const *measuredStarList) const {
36  /**********************************************************************/
37  /* @note the math in this method and accumulateStatImageList() must be kept consistent,
38  * in terms of +/- convention, definition of model, etc. */
39  /**********************************************************************/
40 
41  /* this routine works in two different ways: either providing the
42  Ccd, of providing the MeasuredStarList. In the latter case, the
43  Ccd should match the one(s) in the list. */
44  if (measuredStarList) assert(&(measuredStarList->front()->getCcdImage()) == &ccdImage);
45 
46  auto const &mapping = _photometryModel->getMapping(ccdImage);
47 
48  unsigned nparModel = (_fittingModel) ? mapping.getNpar() : 0;
49  unsigned nparFlux = (_fittingFluxes) ? 1 : 0;
50  unsigned nparTotal = nparModel + nparFlux;
51  std::vector<unsigned> indices(nparModel, -1);
52 
53  Eigen::VectorXd H(nparTotal); // derivative matrix
54  // current position in the Jacobian
55  unsigned kTriplets = tripletList.getNextFreeIndex();
56  const MeasuredStarList &catalog = (measuredStarList) ? *measuredStarList : ccdImage.getCatalogForFit();
57 
58  for (auto const &measuredStar : catalog) {
59  if (!measuredStar->isValid()) continue;
60 
61 // tweak the measurement errors
62 #ifdef FUTURE
63  TweakPhotomMeasurementErrors(inPos, *measuredStar, _fluxError);
64 #endif
65  double fluxErr = measuredStar->getFluxErr();
66  double photomFactor = _photometryModel->photomFactor(ccdImage, *measuredStar);
67  H.setZero(); // we cannot be sure that all entries will be overwritten.
68 
69  double residual = computeMeasurementResidual(photomFactor, *measuredStar);
70 
71  double inverseSigma = 1.0 / (photomFactor * fluxErr);
72  double W = sqr(inverseSigma);
73 
74  if (_fittingModel) {
75  _photometryModel->getMappingIndices(ccdImage, indices);
76  _photometryModel->computeParameterDerivatives(*measuredStar, ccdImage, H);
77  for (unsigned k = 0; k < indices.size(); k++) {
78  unsigned l = indices[k];
79  tripletList.addTriplet(l, kTriplets, H[k] * inverseSigma);
80  grad[l] += H[k] * W * residual;
81  }
82  }
83  if (_fittingFluxes) {
84  unsigned index = measuredStar->getFittedStar()->getIndexInMatrix();
85  // Note: H = dR/dFittedStarFlux == -1
86  tripletList.addTriplet(index, kTriplets, -1.0 * inverseSigma);
87  grad[index] += -1.0 * W * residual;
88  }
89  kTriplets += 1; // each measurement contributes 1 column in the Jacobian
90  }
91 
92  tripletList.setNextFreeIndex(kTriplets);
93 }
94 
95 void PhotometryFit::leastSquareDerivativesReference(FittedStarList const &fittedStarList,
96  TripletList &tripletList, Eigen::VectorXd &grad) const {
97  /**********************************************************************/
100  /**********************************************************************/
101 
102  // Derivatives of terms involving fitted and refstars only contribute if we are fitting fluxes.
103  if (!_fittingFluxes) return;
104  // Can't compute anything if there are no refStars.
105  if (_associations->refStarList.size() == 0) return;
106 
107  unsigned kTriplets = tripletList.getNextFreeIndex();
108 
109  for (auto const &fittedStar : fittedStarList) {
110  auto refStar = fittedStar->getRefStar();
111  if (refStar == nullptr) continue; // no contribution if no associated refstar
112  // TODO: Can we actually work with multiple filters at a time in this scheme?
113  // TODO: I feel we might only be able to fit one filter at a time, because these terms are
114  // independent of the ccdImage, so we don't have a specific filter.
115  // filter = ccdImage.getFilter();
116 
117  // W == inverseSigma^2
118  double inverseSigma = 1.0 / refStar->getFluxErr();
119  // Residual is fittedStar.flux - refStar.flux for consistency with measurement terms.
120  double residual = fittedStar->getFlux() - refStar->getFlux();
121 
122  unsigned index = fittedStar->getIndexInMatrix();
123  // Note: H = dR/dFittedStar == 1
124  tripletList.addTriplet(index, kTriplets, 1.0 * inverseSigma);
125  grad(index) += 1.0 * sqr(inverseSigma) * residual;
126  kTriplets += 1;
127  }
128  tripletList.setNextFreeIndex(kTriplets);
129 }
130 
131 void PhotometryFit::accumulateStatImageList(CcdImageList const &ccdImageList, Chi2Accumulator &accum) const {
132  /**********************************************************************/
135  /**********************************************************************/
136  for (auto const &ccdImage : ccdImageList) {
137  auto &catalog = ccdImage->getCatalogForFit();
138 
139  for (auto const &measuredStar : catalog) {
140  if (!measuredStar->isValid()) continue;
141  double photomFactor = _photometryModel->photomFactor(*ccdImage, *measuredStar);
142  // tweak the measurement errors
143  double sigma = (measuredStar->getFluxErr() * photomFactor);
144 #ifdef FUTURE
145  TweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
146 #endif
147  double residual = computeMeasurementResidual(photomFactor, *measuredStar);
148 
149  double chi2Val = sqr(residual / sigma);
150  accum.addEntry(chi2Val, 1, measuredStar);
151  } // end loop on measurements
152  }
153 }
154 
155 void PhotometryFit::accumulateStatRefStars(Chi2Accumulator &accum) const {
156  /**********************************************************************/
159  /**********************************************************************/
160 
161  FittedStarList &fittedStarList = _associations->fittedStarList;
162  for (auto const &fittedStar : fittedStarList) {
163  auto refStar = fittedStar->getRefStar();
164  if (refStar == nullptr) continue;
165  double chi2 = sqr(((fittedStar->getFlux() - refStar->getFlux()) / refStar->getFluxErr()));
166  accum.addEntry(chi2, 1, fittedStar);
167  }
168 }
169 
171 
173 void PhotometryFit::getIndicesOfMeasuredStar(MeasuredStar const &measuredStar,
174  std::vector<unsigned> &indices) const {
175  indices.clear();
176  if (_fittingModel) {
177  _photometryModel->getMappingIndices(measuredStar.getCcdImage(), indices);
178  }
179  if (_fittingFluxes) {
180  auto fs = measuredStar.getFittedStar();
181  unsigned fsIndex = fs->getIndexInMatrix();
182  indices.push_back(fsIndex);
183  }
184 }
185 
186 void PhotometryFit::assignIndices(std::string const &whatToFit) {
187  _whatToFit = whatToFit;
188  LOGLS_INFO(_log, "assignIndices: now fitting: " << whatToFit);
189  _fittingModel = (_whatToFit.find("Model") != std::string::npos);
190  _fittingFluxes = (_whatToFit.find("Fluxes") != std::string::npos);
191  // When entering here, we assume that whatToFit has already been interpreted.
192 
193  _nParModel = (_fittingModel) ? _photometryModel->assignIndices(whatToFit, 0) : 0;
194  unsigned ipar = _nParModel;
195 
196  if (_fittingFluxes) {
197  for (auto &fittedStar : _associations->fittedStarList) {
198  // the parameter layout here is used also
199  // - when filling the derivatives
200  // - when updating (offsetParams())
201  // - in getIndicesOfMeasuredStar
202  fittedStar->setIndexInMatrix(ipar);
203  ipar += 1;
204  }
205  }
206  _nParFluxes = ipar - _nParModel;
207  _nParTot = ipar;
208  LOGLS_DEBUG(_log,
209  "nParameters total: " << _nParTot << " model: " << _nParModel << " fluxes: " << _nParFluxes);
210 }
211 
212 void PhotometryFit::offsetParams(Eigen::VectorXd const &delta) {
213  if (delta.size() != _nParTot)
214  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
215  "PhotometryFit::offsetParams : the provided vector length is not compatible with "
216  "the current whatToFit setting");
217  if (_fittingModel) _photometryModel->offsetParams(delta);
218 
219  if (_fittingFluxes) {
220  for (auto &fittedStar : _associations->fittedStarList) {
221  // the parameter layout here is used also
222  // - when filling the derivatives
223  // - when assigning indices (assignIndices())
224  unsigned index = fittedStar->getIndexInMatrix();
225  fittedStar->getFlux() -= delta(index);
226  }
227  }
228 }
229 
230 void PhotometryFit::saveResultTuples(std::string const &tupleName) const {
231  std::ofstream tuple(tupleName.c_str());
232  /* If we think the some coordinate on the focal plane is relevant in
233  the ntuple, because thmodel relies on it, then we have to add
234  some function to the model that returns this relevant
235  coordinate. */
236  tuple << "#xccd: coordinate in CCD" << std::endl
237  << "#yccd: " << std::endl
238  << "#mag: rough mag" << std::endl
239  << "#flux : measured flux" << std::endl
240  << "#fluxError : measured flux error" << std::endl
241  << "#fflux : fitted flux" << std::endl
242  << "#phot_factor:" << std::endl
243  << "#jd: Julian date of the measurement" << std::endl
244  << "#color : " << std::endl
245  << "#fsindex: some unique index of the object" << std::endl
246  << "#ra: pos of fitted star" << std::endl
247  << "#dec: pos of fitted star" << std::endl
248  << "#chi2: contribution to Chi2 (1 dof)" << std::endl
249  << "#nm: number of measurements of this FittedStar" << std::endl
250  << "#chip: chip number" << std::endl
251  << "#visit: visit id" << std::endl
252  << "#end" << std::endl;
253  const CcdImageList &ccdImageList = _associations->getCcdImageList();
254  for (auto const &i : ccdImageList) {
255  const CcdImage &im = *i;
256  const MeasuredStarList &cat = im.getCatalogForFit();
257  for (auto const &is : cat) {
258  const MeasuredStar &ms = *is;
259  if (!ms.isValid()) continue;
260  double sigma = ms.getFluxErr();
261 #ifdef FUTURE
262  tweakPhotomMeasurementErrors(inPos, ms, _fluxError);
263 #endif
264  double photomFactor = _photometryModel->photomFactor(im, ms);
265  double jd = im.getMjd();
266  auto fs = ms.getFittedStar();
267  double residual = ms.getFlux() - photomFactor * fs->getFlux();
268  double chi2Val = sqr(residual / sigma);
269  tuple << ms.x << ' ' << ms.y << ' ' << fs->getMag() << ' ' << ms.getFlux() << ' '
270  << ms.getFluxErr() << ' ' << fs->getFlux() << ' ' << photomFactor << ' ' << jd << ' '
271  << fs->color << ' ' << fs->getIndexInMatrix() << ' ' << fs->x << ' ' << fs->y << ' '
272  << chi2Val << ' ' << fs->getMeasurementCount() << ' ' << im.getCcdId() << ' '
273  << im.getVisit() << std::endl;
274  } // loop on measurements in image
275  } // loop on images
276 }
277 } // namespace jointcal
278 } // namespace lsst
int getCcdId() const
returns ccd ID
Definition: CcdImage.h:132
VisitIdType getVisit() const
returns visit ID
Definition: CcdImage.h:135
double getMjd() const
Julian Date.
Definition: CcdImage.h:141
double getFluxErr() const
Definition: BaseStar.h:75
A list of MeasuredStar. They are usually filled in Associations::AddImage.
Definition: MeasuredStar.h:95
double x
coordinate
Definition: Point.h:18
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the requested quantities.
Class for a simple mapping implementing a generic Gtransfo.
Definition: Associations.h:24
void saveResultTuples(std::string const &tupleName) const override
Save the full chi2 term per star that was used in the minimization, for debugging.
objects measured on actual images.
Definition: MeasuredStar.h:18
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:107
double getFlux() const
Definition: BaseStar.h:71
void assignIndices(std::string const &whatToFit) override
Set parameters to fit and assign indices in the big matrix.
std::shared_ptr< const FittedStar > getFittedStar() const
Definition: MeasuredStar.h:80
MeasuredStarList const & getCatalogForFit() const
Gets the catalog to be used for fitting, which may have been cleaned-up.
Definition: CcdImage.h:98
std::list< std::shared_ptr< CcdImage > > CcdImageList
Definition: CcdImage.h:22
Handler of an actual image from a single CCD.
Definition: CcdImage.h:31
bool isValid() const
Fits may use that to discard outliers.
Definition: MeasuredStar.h:87