lsst.meas.algorithms  13.0-24-g22030a45+6
SpatialModelPsf.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 
3 /*
4  * LSST Data Management System
5  * Copyright 2008, 2009, 2010 LSST Corporation.
6  *
7  * This product includes software developed by the
8  * LSST Project (http://www.lsst.org/).
9  *
10  * This program is free software: you can redistribute it and/or modify
11  * it under the terms of the GNU General Public License as published by
12  * the Free Software Foundation, either version 3 of the License, or
13  * (at your option) any later version.
14  *
15  * This program is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the LSST License Statement and
21  * the GNU General Public License along with this program. If not,
22  * see <http://www.lsstcorp.org/LegalNotices/>.
23  */
24 
32 #include <numeric>
33 
34 #if !defined(DOXYGEN)
35 # include "Minuit2/FCNBase.h"
36 # include "Minuit2/FunctionMinimum.h"
37 # include "Minuit2/MnMigrad.h"
38 # include "Minuit2/MnMinos.h"
39 # include "Minuit2/MnPrint.h"
40 #endif
41 
42 #include "Eigen/Core"
43 #include "Eigen/Cholesky"
44 #include "Eigen/SVD"
45 
46 #include "lsst/afw/detection/Footprint.h"
47 #include "lsst/afw/math/SpatialCell.h"
48 #include "lsst/afw/math/FunctionLibrary.h"
49 #include "lsst/afw/geom/Point.h"
50 #include "lsst/afw/geom/Box.h"
54 
56 namespace afwGeom = lsst::afw::geom;
57 namespace afwImage = lsst::afw::image;
58 namespace afwMath = lsst::afw::math;
59 
60 namespace lsst {
61 namespace meas {
62 namespace algorithms {
63 
64 namespace {
65 
66 int const WARP_BUFFER(1); // Buffer (border) around kernel to prevent warp issues
67 std::string const WARP_ALGORITHM("lanczos5"); // Name of warping algorithm to use
68 
69 
70 // A class to pass around to all our PsfCandidates which builds the PcaImageSet
71 template<typename PixelT>
72 class SetPcaImageVisitor : public afwMath::CandidateVisitor {
73  typedef afwImage::Image<PixelT> ImageT;
74  typedef afwImage::MaskedImage<PixelT> MaskedImageT;
75  typedef afwImage::Exposure<PixelT> ExposureT;
76 public:
77  explicit SetPcaImageVisitor(
78  PsfImagePca<MaskedImageT> *imagePca, // Set of Images to initialise
79  unsigned int const mask=0x0 // Ignore pixels with any of these bits set
80  ) :
81  afwMath::CandidateVisitor(),
82  _imagePca(imagePca)
83  {
84  ;
85  }
86 
87  // Called by SpatialCellSet::visitCandidates for each Candidate
88  void processCandidate(afwMath::SpatialCellCandidate *candidate) {
89  PsfCandidate<PixelT> *imCandidate = dynamic_cast<PsfCandidate<PixelT> *>(candidate);
90  if (imCandidate == NULL) {
91  throw LSST_EXCEPT(lsst::pex::exceptions::LogicError,
92  "Failed to cast SpatialCellCandidate to PsfCandidate");
93  }
94 
95  try {
96  std::shared_ptr<MaskedImageT> im = imCandidate->getOffsetImage(WARP_ALGORITHM,
97  WARP_BUFFER);
98 
99 
100  //static int count = 0;
101  //im->writeFits(str(boost::format("cand%03d.fits") % count));
102  //count += 1;
103 
104  afwMath::StatisticsControl sctrl;
105  sctrl.setNanSafe(false);
106 
107  if (!std::isfinite(afwMath::makeStatistics(*im->getImage(),
108  afwMath::MAX, sctrl).getValue())) {
109  throw LSST_EXCEPT(lsst::pex::exceptions::RuntimeError,
110  str(boost::format("Image at %d, %d contains NaN")
111  % imCandidate->getXCenter() % imCandidate->getYCenter()));
112 
113  }
114  if (!std::isfinite(afwMath::makeStatistics(*im->getVariance(),
115  afwMath::MAX, sctrl).getValue())) {
116  throw LSST_EXCEPT(lsst::pex::exceptions::RuntimeError,
117  str(boost::format("Variance of Image at %d, %d contains NaN")
118  % imCandidate->getXCenter() % imCandidate->getYCenter()));
119  }
120 
121  _imagePca->addImage(im, imCandidate->getSource()->getPsfFlux());
122  } catch(lsst::pex::exceptions::LengthError &) {
123  return;
124  }
125  }
126 private:
127  PsfImagePca<MaskedImageT> *_imagePca; // the ImagePca we're building
128 };
129 
130 /************************************************************************************************************/
132 template<typename PixelT>
133 class countVisitor : public afwMath::CandidateVisitor {
134  typedef afwImage::MaskedImage<PixelT> MaskedImage;
135  typedef afwImage::Exposure<PixelT> Exposure;
136 public:
137  explicit countVisitor() : afwMath::CandidateVisitor(), _n(0) {}
138 
139  void reset() {
140  _n = 0;
141  }
142 
143  // Called by SpatialCellSet::visitCandidates for each Candidate
144  void processCandidate(afwMath::SpatialCellCandidate *candidate) {
145  PsfCandidate<PixelT> *imCandidate = dynamic_cast<PsfCandidate<PixelT> *>(candidate);
146  if (imCandidate == NULL) {
147  throw LSST_EXCEPT(lsst::pex::exceptions::LogicError,
148  "Failed to cast SpatialCellCandidate to PsfCandidate");
149  }
150 
151  try {
152  imCandidate->getMaskedImage();
153  } catch(lsst::pex::exceptions::LengthError &) {
154  return;
155  }
156 
157  ++_n;
158  }
159 
160  // Return the number
161  double getN() const { return _n; }
162 
163 private:
164  int mutable _n; // the desired number
165 };
166 
167 
172 template<typename ImageT>
173 std::vector<std::shared_ptr<ImageT>> offsetKernel(
174  afwMath::LinearCombinationKernel const& kernel,
175  float dx, float dy
176  )
177 {
178  afwMath::KernelList kernels = kernel.getKernelList(); // The Kernels that kernel adds together
179  unsigned int const nKernel = kernels.size(); // Number of kernel components
180  std::vector<std::shared_ptr<ImageT>> kernelImages(nKernel); // Images of each Kernel in kernels
181  if (nKernel == 0) {
182  throw LSST_EXCEPT(lsst::pex::exceptions::LengthError,
183  "Kernel has no components");
184  }
185 
186  ImageT scratch(kernel.getDimensions()); // Buffered scratch space
187  for (unsigned int i = 0; i != nKernel; ++i) {
188  kernels[i]->computeImage(scratch, false);
189  kernelImages[i] = afwMath::offsetImage(scratch, dx, dy, WARP_ALGORITHM, WARP_BUFFER);
190  }
191 
192  return kernelImages;
193 }
194 
195 } // Anonymous namespace
196 
197 /************************************************************************************************************/
205 template<typename PixelT>
206 std::pair<std::shared_ptr<afwMath::LinearCombinationKernel>, std::vector<double> > createKernelFromPsfCandidates(
207  afwMath::SpatialCellSet const& psfCells,
208  lsst::afw::geom::Extent2I const& dims,
209  lsst::afw::geom::Point2I const& xy0,
210  int const nEigenComponents,
211  int const spatialOrder,
212  int const ksize,
213  int const nStarPerCell,
214  bool const constantWeight,
215  int const border
216  )
217 {
218  typedef typename afwImage::Image<PixelT> ImageT;
219  typedef typename afwImage::MaskedImage<PixelT> MaskedImageT;
220  typedef typename afwImage::Exposure<PixelT> ExposureT;
221 
222  //
223  // Set the sizes for PsfCandidates made from either Images or MaskedImages
224  //
225  //lsst::meas::algorithms::PsfCandidate<ImageT>::setWidth(ksize);
226  //lsst::meas::algorithms::PsfCandidate<ImageT>::setHeight(ksize);
227  //lsst::meas::algorithms::PsfCandidate<MaskedImageT>::setWidth(ksize);
228  //lsst::meas::algorithms::PsfCandidate<MaskedImageT>::setHeight(ksize);
231 
232 
233  PsfImagePca<MaskedImageT> imagePca(constantWeight, border); // Here's the set of images we'll analyze
234 
235  {
236  SetPcaImageVisitor<PixelT> importStarVisitor(&imagePca);
237  bool const ignoreExceptions = true;
238  psfCells.visitCandidates(&importStarVisitor, nStarPerCell, ignoreExceptions);
239  }
240 
241  //
242  // Do a PCA decomposition of those PSF candidates.
243  //
244  // We have "gappy" data; in other words we don't want to include any pixels with INTRP set
245  //
246  int niter = 10; // number of iterations of updateBadPixels
247  double deltaLim = 10.0; // acceptable value of delta, the max change due to updateBadPixels
248  lsst::afw::image::MaskPixel const BAD = afwImage::Mask<>::getPlaneBitMask("BAD");
249  lsst::afw::image::MaskPixel const CR = afwImage::Mask<>::getPlaneBitMask("CR");
250  lsst::afw::image::MaskPixel const INTRP = afwImage::Mask<>::getPlaneBitMask("INTRP");
251 
252  for (int i = 0; i != niter; ++i) {
253  int const ncomp = (i == 0) ? 0 :
254  ((nEigenComponents == 0) ? imagePca.getEigenImages().size() : nEigenComponents);
255  double delta = imagePca.updateBadPixels(BAD | CR | INTRP, ncomp);
256  if (i > 0 && delta < deltaLim) {
257  break;
258  }
259 
260  imagePca.analyze();
261  }
262 
263  std::vector<std::shared_ptr<MaskedImageT>> eigenImages = imagePca.getEigenImages();
264  std::vector<double> eigenValues = imagePca.getEigenValues();
265  int const nEigen = static_cast<int>(eigenValues.size());
266 
267  int const ncomp = (nEigenComponents <= 0 || nEigen < nEigenComponents) ? nEigen : nEigenComponents;
268  //
269  // Set the background level of the components to 0.0 to avoid coupling variable background
270  // levels to the form of the Psf. More precisely, we calculate the mean of an outer "annulus"
271  // of width bkg_border
272  //
273  for (int k = 0; k != ncomp; ++k) {
274  ImageT const& im = *eigenImages[k]->getImage();
275 
276  int bkg_border = 2;
277  if (bkg_border > im.getWidth()) {
278  bkg_border = im.getWidth() / 2;
279  }
280  if (bkg_border > im.getHeight()) {
281  bkg_border = im.getHeight() / 2;
282  }
283 
284  double sum = 0;
285  // Bottom and Top borders
286  for (int i = 0; i != bkg_border; ++i) {
287  typename ImageT::const_x_iterator
288  ptrB = im.row_begin(i), ptrT = im.row_begin(im.getHeight() - i - 1);
289  for (int j = 0; j != im.getWidth(); ++j, ++ptrB, ++ptrT) {
290  sum += *ptrB + *ptrT;
291  }
292  }
293  for (int i = bkg_border; i < im.getHeight() - bkg_border; ++i) {
294  // Left and Right borders
295  typename ImageT::const_x_iterator
296  ptrL = im.row_begin(i), ptrR = im.row_begin(i) + im.getWidth() - bkg_border;
297  for (int j = 0; j != bkg_border; ++j, ++ptrL, ++ptrR) {
298  sum += *ptrL + *ptrR;
299  }
300  }
301  sum /= 2*(bkg_border*im.getWidth() + bkg_border*(im.getHeight() - 2*bkg_border));
302 
303  *eigenImages[k] -= sum;
304  }
305  //
306  // Now build our LinearCombinationKernel; build the lists of basis functions
307  // and spatial variation, then assemble the Kernel
308  //
309  afwMath::KernelList kernelList;
310  std::vector<afwMath::Kernel::SpatialFunctionPtr> spatialFunctionList;
311  afwGeom::Box2D const range = afwGeom::Box2D(afwGeom::Point2D(xy0), afwGeom::Extent2D(dims));
312 
313  for (int i = 0; i != ncomp; ++i) {
314  {
315  // Enforce unit sum for kernel by construction
316  // Zeroth component has unit sum
317  // Other components have zero sum by normalising and then subtracting the zeroth component
318  ImageT& image = *eigenImages[i]->getImage();
319  double sum = std::accumulate(image.begin(true), image.end(true), 0.0);
320  if (i == 0) {
321  image /= sum;
322  } else {
323  for (typename ImageT::fast_iterator ptr0 = eigenImages[0]->getImage()->begin(true),
324  ptr1 = image.begin(true), end = image.end(true); ptr1 != end; ++ptr0, ++ptr1) {
325  *ptr1 = *ptr1 / sum - *ptr0;
326  }
327  }
328  }
329 
330  kernelList.push_back(std::shared_ptr<afwMath::Kernel>(new afwMath::FixedKernel(
331  afwImage::Image<afwMath::Kernel::Pixel>(*eigenImages[i]->getImage(),true)
332  )));
333 
334  afwMath::Kernel::SpatialFunctionPtr
335 // spatialFunction(new afwMath::PolynomialFunction2<double>(spatialOrder));
336  spatialFunction(new afwMath::Chebyshev1Function2<double>(spatialOrder, range));
337  spatialFunction->setParameter(0, 1.0); // the constant term; all others are 0
338  spatialFunctionList.push_back(spatialFunction);
339  }
340 
341  std::shared_ptr<afwMath::LinearCombinationKernel>
342  psf(new afwMath::LinearCombinationKernel(kernelList, spatialFunctionList));
343 
344  return std::make_pair(psf, eigenValues);
345 }
346 
347 /************************************************************************************************************/
351 template<typename PixelT>
352 int countPsfCandidates(afwMath::SpatialCellSet const& psfCells,
353  int const nStarPerCell)
354 {
355  countVisitor<PixelT> counter;
356  psfCells.visitCandidates(&counter, nStarPerCell);
357 
358  return counter.getN();
359 }
360 
361 /************************************************************************************************************/
362 namespace {
368 template<typename ModelImageT, typename DataImageT>
369 std::pair<double, double>
370 fitKernel(ModelImageT const& mImage, // The model image at this point
371  DataImageT const& data, // the data to fit
372  double lambda = 0.0, // floor for variance is lambda*data
373  bool detected = true, // only fit DETECTED pixels?
374  int const id=-1 // ID for this object; useful in debugging
375  ) {
376  assert(data.getDimensions() == mImage.getDimensions());
377  assert(id == id);
378  int const DETECTED = afwImage::Mask<>::getPlaneBitMask("DETECTED");
379  int const BAD = afwImage::Mask<>::getPlaneBitMask("CR") | afwImage::Mask<>::getPlaneBitMask("BAD");
380 
381  double sumMM = 0.0, sumMD = 0.0, sumDD = 0.0; // sums of model*model/variance etc.
382  int npix = 0; // number of pixels used to evaluate chi^2
383  for (int y = 0; y != data.getHeight(); ++y) {
384  typename ModelImageT::x_iterator mptr = mImage.row_begin(y);
385  for (typename DataImageT::x_iterator ptr = data.row_begin(y), end = data.row_end(y);
386  ptr != end; ++ptr, ++mptr) {
387  double const m = (*mptr)[0]; // value of model
388  double const d = ptr.image(); // value of data
389  double const var = ptr.variance() + lambda*d; // data's variance
390  if (detected && !(ptr.mask() & DETECTED)) {
391  continue;
392  }
393  if (ptr.mask() & BAD) {
394  continue;
395  }
396  if (var != 0.0) { // assume variance == 0 => infinity XXX
397  double const iVar = 1.0/var;
398  npix++;
399  sumMM += m*m*iVar;
400  sumMD += m*d*iVar;
401  sumDD += d*d*iVar;
402  }
403  }
404  }
405 
406  if (npix == 0) {
407  throw LSST_EXCEPT(lsst::pex::exceptions::RangeError, "No good pixels");
408  }
409  if (sumMM == 0.0) {
410  throw LSST_EXCEPT(lsst::pex::exceptions::RangeError, "sum(data*data)/var == 0");
411  }
412 
413  double const amp = sumMD/sumMM; // estimate of amplitude of model at this point
414  double const chi2 = (sumDD - 2*amp*sumMD + amp*amp*sumMM)/(npix - 1);
415 
416 #if 0
417  bool show = false; // Display the centre of the image; set from gdb
418 
419  if (show) {
420  show = true; // you can jump here in gdb to set show if direct attempts fail
421  int y = data.getHeight()/2;
422  int x = data.getWidth()/2;
423  int hsize = 2;
424  printf("\ndata ");
425  for (int ii = -hsize; ii <= hsize; ++ii) {
426  for (int jj = -hsize; jj <= hsize; ++jj) {
427  printf("%7.1f ", data.at(x + jj, y - ii).image());
428  }
429  printf(" model ");
430  for (int jj = -hsize; jj <= hsize; ++jj) {
431  printf("%7.1f ", amp*(*(mImage.at(x + jj, y - ii)))[0]);
432  }
433  printf("\n ");
434  }
435  printf("%g %.1f\n", amp, chi2);
436  }
437 #endif
438 
439  return std::make_pair(chi2, amp);
440 }
441 }
442 
443 /************************************************************************************************************/
444 /*
445  * Fit for the spatial variation of the PSF parameters over the field
446  */
448 template<typename PixelT>
449 class evalChi2Visitor : public afwMath::CandidateVisitor {
451  typedef afwImage::MaskedImage<PixelT> MaskedImage;
452  typedef afwImage::Exposure<PixelT> Exposure;
453 
455 public:
456  explicit evalChi2Visitor(afwMath::Kernel const& kernel,
457  double lambda
458  ) :
459  afwMath::CandidateVisitor(),
460  _chi2(0.0), _kernel(kernel), _lambda(lambda),
461  _kImage(std::shared_ptr<KImage>(new KImage(kernel.getDimensions()))) {
462  }
463 
464  void reset() {
465  _chi2 = 0.0;
466  }
467 
468  // Called by SpatialCellSet::visitCandidates for each Candidate
469  void processCandidate(afwMath::SpatialCellCandidate *candidate) {
470  PsfCandidate<PixelT> *imCandidate = dynamic_cast<PsfCandidate<PixelT> *>(candidate);
471  if (imCandidate == NULL) {
472  throw LSST_EXCEPT(lsst::pex::exceptions::LogicError,
473  "Failed to cast SpatialCellCandidate to PsfCandidate");
474  }
475 
476  double const xcen = imCandidate->getSource()->getX();
477  double const ycen = imCandidate->getSource()->getY();
478 
479  _kernel.computeImage(*_kImage, true, xcen, ycen);
480  std::shared_ptr<MaskedImage const> data;
481  try {
482  data = imCandidate->getOffsetImage(WARP_ALGORITHM, WARP_BUFFER);
483  } catch(lsst::pex::exceptions::LengthError &) {
484  return;
485  }
486 
487  try {
488  std::pair<double, double> result = fitKernel(*_kImage, *data, _lambda, false,
489  imCandidate->getSource()->getId());
490 
491  double dchi2 = result.first; // chi^2 from this object
492  double const amp = result.second; // estimate of amplitude of model at this point
493 
494  imCandidate->setChi2(dchi2);
495  imCandidate->setAmplitude(amp);
496 
497  _chi2 += dchi2;
498  } catch(lsst::pex::exceptions::RangeError &e) {
499  imCandidate->setStatus(afwMath::SpatialCellCandidate::BAD);
500  imCandidate->setChi2(std::numeric_limits<double>::quiet_NaN());
501  imCandidate->setAmplitude(std::numeric_limits<double>::quiet_NaN());
502  }
503  }
504 
505  // Return the computed chi^2
506  double getValue() const { return _chi2; }
507 
508 private:
509  double mutable _chi2; // the desired chi^2
510  afwMath::Kernel const& _kernel; // the kernel
511  double _lambda; // floor for variance is _lambda*data
512  std::shared_ptr<KImage> mutable _kImage; // The Kernel at this point; a scratch copy
513 };
514 
515 /********************************************************************************************************/
519 // Set the Kernel's spatial parameters from a vector
520 void setSpatialParameters(afwMath::Kernel *kernel,
521  std::vector<double> const& coeffs
522  )
523 {
524  int const nComponents = kernel->getNKernelParameters();
525  int const nSpatialParams = kernel->getNSpatialParameters();
526 
527  assert (nComponents*nSpatialParams == static_cast<long>(coeffs.size()));
528 
529  std::vector<std::vector<double> > kCoeffs; // coefficients rearranged for Kernel
530  kCoeffs.reserve(nComponents);
531  for (int i = 0; i != nComponents; ++i) {
532  kCoeffs.push_back(std::vector<double>(nSpatialParams));
533  std::copy(coeffs.begin() + i*nSpatialParams,
534  coeffs.begin() + (i + 1)*nSpatialParams, kCoeffs[i].begin());
535  }
536 
537  kernel->setSpatialParameters(kCoeffs);
538 }
539 
543 // Set the Kernel's spatial parameters from an Eigen::VectorXd
544 void setSpatialParameters(afwMath::Kernel *kernel,
545  Eigen::VectorXd const& vec
546  )
547 {
548  int const nComponents = kernel->getNKernelParameters();
549  int const nSpatialParams = kernel->getNSpatialParameters();
550 
551  assert (nComponents*nSpatialParams == vec.size());
552 
553  std::vector<std::vector<double> > kCoeffs; // coefficients rearranged for Kernel
554  kCoeffs.reserve(nComponents);
555  for (int i = 0; i != nComponents; ++i) {
556  std::vector<double> spatialCoeffs(nSpatialParams);
557  for (int j = 0; j != nSpatialParams; ++j) {
558  spatialCoeffs[j] = vec[i*nSpatialParams + j];
559  }
560  kCoeffs.push_back(spatialCoeffs);
561  }
562 
563  kernel->setSpatialParameters(kCoeffs);
564 }
565 
566 //
567 // The object that minuit minimises
568 //
569 template<typename PixelT>
570 class MinimizeChi2 : public ROOT::Minuit2::FCNBase {
571 public:
572  explicit MinimizeChi2(evalChi2Visitor<PixelT> & chi2Visitor,
573  afwMath::Kernel *kernel,
574  afwMath::SpatialCellSet const& psfCells,
575  int nStarPerCell,
576  int nComponents,
577  int nSpatialParams
578  ) : _errorDef(1.0),
579  _chi2Visitor(chi2Visitor),
580  _kernel(kernel),
581  _psfCells(psfCells),
582  _nStarPerCell(nStarPerCell),
583  _nComponents(nComponents),
584  _nSpatialParams(nSpatialParams) {}
585 
592  double Up() const { return _errorDef; }
593 
594  // Evaluate our cost function (in this case chi^2)
595  double operator()(const std::vector<double>& coeffs) const {
596  setSpatialParameters(_kernel, coeffs);
597 
598  _psfCells.visitCandidates(&_chi2Visitor, _nStarPerCell);
599 
600  return _chi2Visitor.getValue();
601  }
602 
603  void setErrorDef(double def) { _errorDef = def; }
604 private:
605  double _errorDef; // how much cost function has changed at the +- 1 error points
606 
607  evalChi2Visitor<PixelT>& _chi2Visitor;
608  afwMath::Kernel *_kernel;
609  afwMath::SpatialCellSet const& _psfCells;
610  int _nStarPerCell;
611  int _nComponents;
612  int _nSpatialParams;
613 };
614 
615 /************************************************************************************************************/
619 template<typename PixelT>
620 std::pair<bool, double>
622  afwMath::Kernel *kernel,
623  afwMath::SpatialCellSet const& psfCells,
624  int const nStarPerCell,
625  double const tolerance,
626  double const lambda
627  ) {
628  typedef typename afwImage::Image<PixelT> Image;
629 
630  int const nComponents = kernel->getNKernelParameters();
631  int const nSpatialParams = kernel->getNSpatialParameters();
632  //
633  // visitor that evaluates the chi^2 of the current fit
634  //
635  evalChi2Visitor<PixelT> getChi2(*kernel, lambda);
636  //
637  // We have to unpack the Kernel coefficients into a linear array, coeffs
638  //
639  std::vector<double> coeffs; // The coefficients we want to fit
640  coeffs.assign(nComponents*nSpatialParams, 0.0);
641 
642  std::vector<double> stepSize; // step sizes
643  stepSize.assign(nComponents*nSpatialParams, 100);
644  //
645  // Translate that into minuit's language
646  //
647  ROOT::Minuit2::MnUserParameters fitPar;
648  std::vector<std::string> paramNames;
649  paramNames.reserve(nComponents*nSpatialParams);
650 
651  for (int i = 0, c = 0; c != nComponents; ++c) {
652  coeffs[i] = 1; // the constant part of each spatial order
653  for (int s = 0; s != nSpatialParams; ++s, ++i) {
654  paramNames.push_back((boost::format("C%d:%d") % c % s).str());
655  fitPar.Add(paramNames[i].c_str(), coeffs[i], stepSize[i]);
656  }
657  }
658  fitPar.Fix("C0:0");
659  //
660  // Create the minuit object that knows how to minimise our functor
661  //
662  MinimizeChi2<PixelT> minimizerFunc(getChi2, kernel, psfCells, nStarPerCell, nComponents, nSpatialParams);
663 
664  double const errorDef = 1.0; // use +- 1sigma errors
665  minimizerFunc.setErrorDef(errorDef);
666  //
667  // tell minuit about it
668  //
669  ROOT::Minuit2::MnMigrad migrad(minimizerFunc, fitPar);
670  //
671  // And let it loose
672  //
673  int maxFnCalls = 0; // i.e. unlimited
674  ROOT::Minuit2::FunctionMinimum min =
675  migrad(maxFnCalls, tolerance/(1e-4*errorDef)); // minuit uses 0.1*1e-3*tolerance*errorDef
676 
677  float minChi2 = min.Fval();
678  bool const isValid = min.IsValid() && std::isfinite(minChi2);
679 
680  if (true || isValid) { // calculate coeffs even in minuit is unhappy
681  for (int i = 0; i != nComponents*nSpatialParams; ++i) {
682  coeffs[i] = min.UserState().Value(i);
683  }
684 
685  setSpatialParameters(kernel, coeffs);
686  }
687 
688 #if 0 // Estimate errors; we don't really need this
689  ROOT::Minuit2::MnMinos minos(minimizerFunc, min);
690  for (int i = 0, c = 0; c != nComponents; ++c) {
691  for (int s = 0; s != nSpatialParams; ++s, ++i) {
692  char const *name = paramNames[i].c_str();
693  printf("%s %g", name, min.UserState().Value(name));
694  if (isValid && !fitPar.Parameter(fitPar.Index(name)).IsFixed()) {
695  printf(" (%g+%g)\n", minos(i).first, minos(i).second);
696  }
697  printf("\n");
698  }
699  }
700 #endif
701  //
702  // One time more through the Candidates setting their chi^2 values. We'll
703  // do all the candidates this time, not just the first nStarPerCell
704  //
705  psfCells.visitAllCandidates(&getChi2, true);
706 
707  return std::make_pair(isValid, minChi2);
708 }
709 
710 /************************************************************************************************************/
714 namespace {
751 template<typename PixelT>
752 class FillABVisitor : public afwMath::CandidateVisitor {
753  typedef afwImage::Image<PixelT> Image;
754  typedef afwImage::MaskedImage<PixelT> MaskedImage;
755  typedef afwImage::Exposure<PixelT> Exposure;
756 
758 public:
759  explicit FillABVisitor(afwMath::LinearCombinationKernel const& kernel, // the Kernel we're fitting
760  double tau2=0.0 // floor to the per-candidate variance
761  ) :
762  afwMath::CandidateVisitor(),
763  _kernel(kernel),
764  _tau2(tau2),
765  _nSpatialParams(_kernel.getNSpatialParameters()),
766  _nComponents(_kernel.getNKernelParameters()),
767  _basisImgs(),
768  _A((_nComponents-1)*_nSpatialParams, (_nComponents-1)*_nSpatialParams),
769  _b((_nComponents-1)*_nSpatialParams),
770  _basisDotBasis(_nComponents, _nComponents)
771  {
772  _basisImgs.resize(_nComponents);
773 
774  _A.setZero();
775  _b.setZero();
776  //
777  // Get all the Kernel's components as Images
778  //
779  afwMath::KernelList const& kernels = _kernel.getKernelList(); // Kernel's components
780  for (int i = 0; i != _nComponents; ++i) {
781  _basisImgs[i] = std::shared_ptr<KImage>(new KImage(kernels[i]->getDimensions()));
782  kernels[i]->computeImage(*_basisImgs[i], false);
783  }
784 
785  //
786  // Calculate the inner products of the Kernel components once and for all
787  //
788  for (int i = 1; i != _nComponents; ++i) { // Don't need 0th component
789  for (int j = i; j != _nComponents; ++j) {
790  _basisDotBasis(i, j) = _basisDotBasis(j, i) =
791  afwImage::innerProduct(*_basisImgs[i], *_basisImgs[j],
793  }
794  }
795  }
796 
797  void reset() {}
798 
799  // Called by SpatialCellSet::visitCandidates for each Candidate
800  void processCandidate(afwMath::SpatialCellCandidate *candidate) {
801  PsfCandidate<PixelT> *imCandidate = dynamic_cast<PsfCandidate<PixelT> *>(candidate);
802  if (imCandidate == NULL) {
803  throw LSST_EXCEPT(lsst::pex::exceptions::LogicError,
804  "Failed to cast SpatialCellCandidate to PsfCandidate");
805  }
806 
807  CONST_PTR(MaskedImage) data;
808  try {
809  data = imCandidate->getMaskedImage(_kernel.getWidth(), _kernel.getHeight());
810  } catch(lsst::pex::exceptions::LengthError &) {
811  return;
812  }
813  double const xcen = imCandidate->getXCenter();
814  double const ycen = imCandidate->getYCenter();
815  double const dx = afwImage::positionToIndex(xcen, true).second;
816  double const dy = afwImage::positionToIndex(ycen, true).second;
817 
818 #if 0
819  double const amp = imCandidate->getAmplitude();
820 #else
821  /*
822  * Estimate the amplitude based on the current basis functions.
823  *
824  * N.b. you have to be a little careful here. Consider a PSF that is phi == (N0 + b*y*N1)/(1 + b*y)
825  * where the amplitude of N0 and N1 is 1.0, so a star has profile I = A*(N0 + b*y*N1)/(1 + b*y)
826  *
827  * If we set the amplitude to be A = I(0)/phi(0) (i.e. the central value of the data and best-fit phi)
828  * then the coefficient of N0 becomes 1/(1 + b*y) which makes the model non-linear in y.
829  */
830  std::pair<std::shared_ptr<afwMath::Kernel>, std::pair<double, double> > ret =
831  fitKernelToImage(_kernel, *data, afwGeom::Point2D(xcen, ycen));
832  double const amp = ret.second.first;
833 #endif
834 
835  double const var = imCandidate->getVar();
836  double const ivar = 1/(var + _tau2); // Allow for floor on variance
837 
838  // Spatial params of all the components
839  std::vector<std::vector<double> > params(_nComponents);
840  for (int ic = 1; ic != _nComponents; ++ic) { // Don't need params[0]
841  params[ic] = _kernel.getSpatialFunction(ic)->getDFuncDParameters(xcen, ycen);
842  }
843 
844  std::vector<std::shared_ptr<KImage>> basisImages = offsetKernel<KImage>(_kernel, dx, dy);
845 
846  // Prepare values for basis dot data
847  // Scale data and subtract 0th component as part of unit kernel sum construction
848  std::shared_ptr<Image> dataImage(new Image(*data->getImage(), true));
849  typename KImage::fast_iterator bPtr = basisImages[0]->begin(true);
850  for (typename Image::fast_iterator dPtr = dataImage->begin(true), end = dataImage->end(true);
851  dPtr != end; ++dPtr, ++bPtr) {
852  *dPtr = *dPtr / amp - *bPtr;
853  }
854 
855  for (int i = 0, ic = 1; ic != _nComponents; ++ic) { // Don't need 0th component now
856  double const basisDotData = afwImage::innerProduct(*basisImages[ic], *dataImage,
858  for (int is = 0; is != _nSpatialParams; ++is, ++i) {
859  _b(i) += ivar*params[ic][is]*basisDotData;
860 
861  for (int j = i, jc = ic; jc != _nComponents; ++jc) {
862  for (int js = (i == j) ? is : 0; js != _nSpatialParams; ++js, ++j) {
863  _A(i, j) += ivar*params[ic][is]*params[jc][js]*_basisDotBasis(ic, jc);
864  _A(j, i) = _A(i, j); // could do this after _A is fully calculated
865  }
866  }
867  }
868  }
869  }
870 
871  Eigen::MatrixXd const& getA() const { return _A; }
872  Eigen::VectorXd const& getB() const { return _b; }
873 
874 private:
875  afwMath::LinearCombinationKernel const& _kernel; // the kernel
876  double _tau2; // variance floor added in quadrature to true candidate variance
877  int const _nSpatialParams; // number of spatial parameters
878  int const _nComponents; // number of basis functions
879  std::vector<std::shared_ptr<KImage>> _basisImgs; // basis function images from _kernel
880  Eigen::MatrixXd _A; // We'll solve the matrix equation A x = b for the Kernel's coefficients
881  Eigen::VectorXd _b;
882  Eigen::MatrixXd _basisDotBasis; // the inner products of the Kernel components
883 };
884 
885 
887 template<typename PixelT>
888 class setAmplitudeVisitor : public afwMath::CandidateVisitor {
889  typedef afwImage::MaskedImage<PixelT> MaskedImage;
890  typedef afwImage::Exposure<PixelT> Exposure;
891 public:
892  // Called by SpatialCellSet::visitCandidates for each Candidate
893  void processCandidate(afwMath::SpatialCellCandidate *candidate) {
894  PsfCandidate<PixelT> *imCandidate = dynamic_cast<PsfCandidate<PixelT> *>(candidate);
895  if (imCandidate == NULL) {
896  throw LSST_EXCEPT(lsst::pex::exceptions::LogicError,
897  "Failed to cast SpatialCellCandidate to PsfCandidate");
898  }
899  imCandidate->setAmplitude(afwMath::makeStatistics(*imCandidate->getMaskedImage()->getImage(),
900  afwMath::MAX).getValue());
901  }
902 };
903 
904 }
905 
906 
907 template<typename PixelT>
908 std::pair<bool, double>
910  afwMath::Kernel *kernel,
911  afwMath::SpatialCellSet const& psfCells,
912  bool const doNonLinearFit,
913  int const nStarPerCell,
914  double const tolerance,
915  double const lambda
916  )
917 {
918  if (doNonLinearFit) {
919  return fitSpatialKernelFromPsfCandidates<PixelT>(kernel, psfCells, nStarPerCell, tolerance);
920  }
921 
922  double const tau = 0; // softening for errors
923 
924  afwMath::LinearCombinationKernel const* lcKernel =
925  dynamic_cast<afwMath::LinearCombinationKernel const*>(kernel);
926  if (!lcKernel) {
927  throw LSST_EXCEPT(lsst::pex::exceptions::LogicError,
928  "Failed to cast Kernel to LinearCombinationKernel while building spatial PSF model");
929  }
930 #if 1
931  //
932  // Set the initial amplitudes of all our candidates
933  //
934  setAmplitudeVisitor<PixelT> setAmplitude;
935  psfCells.visitAllCandidates(&setAmplitude, true);
936 #endif
937  //
938  // visitor that fills out the A and b matrices (we'll solve A x = b for the coeffs, x)
939  //
940  FillABVisitor<PixelT> getAB(*lcKernel, tau);
941  //
942  // Actually visit all our candidates
943  //
944  psfCells.visitCandidates(&getAB, nStarPerCell, true);
945  //
946  // Extract A and b, and solve Ax = b
947  //
948  Eigen::MatrixXd const& A = getAB.getA();
949  Eigen::VectorXd const& b = getAB.getB();
950  Eigen::VectorXd x0(b.size()); // Solution to matrix problem
951 
952  switch (b.size()) {
953  case 0: // One candidate, no spatial variability
954  break;
955  case 1: // eigen can't/won't handle 1x1 matrices
956  x0(0) = b(0)/A(0, 0);
957  break;
958  default:
959  x0 = A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
960  break;
961  }
962 #if 0
963  std::cout << "A " << A << std::endl;
964  std::cout << "b " << b.transpose() << std::endl;
965  std::cout << "x " << x.transpose() << std::endl;
966 
967  afwImage::Image<double> img(b.size(), b.size());
968  for (int j = 0; j < b.size(); ++j) {
969  for (int i = 0; i < b.size(); ++i) {
970  img(i, j) = A(i, j);
971  }
972  }
973  img.writeFits("a.fits");
974 
975  if (x.cols() >= 6) {
976  for (int i = 0; i != 6; ++i) {
977  double xcen = 25; double ycen = 35 + 35*i;
978  std::cout << "x, y " << xcen << " , " << ycen << " b "
979  << (x[3] + xcen*x[4] + ycen*x[5])/(x[0] + xcen*x[1] + ycen*x[2]) << std::endl;
980  }
981  }
982 #endif
983 
984  // Generate kernel parameters (including 0th component) from matrix solution
985  Eigen::VectorXd x(kernel->getNKernelParameters() * kernel->getNSpatialParameters()); // Kernel parameters
986  x(0) = 1.0;
987  std::fill(x.data() + 1, x.data() + kernel->getNSpatialParameters(), 0.0);
988  std::copy(x0.data(), x0.data() + x0.size(), x.data() + kernel->getNSpatialParameters());
989 
990  setSpatialParameters(kernel, x);
991  //
992  // One time more through the Candidates setting their chi^2 values. We'll
993  // do all the candidates this time, not just the first nStarPerCell
994  //
995  // visitor that evaluates the chi^2 of the current fit
996  //
997  evalChi2Visitor<PixelT> getChi2(*kernel, lambda);
998 
999  psfCells.visitAllCandidates(&getChi2, true);
1000 
1001  return std::make_pair(true, getChi2.getValue());
1002 }
1003 
1004 /************************************************************************************************************/
1008 template<typename MaskedImageT>
1009 double subtractPsf(afwDetection::Psf const& psf,
1010  MaskedImageT *data,
1011  double x,
1012  double y,
1013  double psfFlux
1014  )
1015 {
1016  if (std::isnan(x + y)) {
1017  return std::numeric_limits<double>::quiet_NaN();
1018  }
1019 
1020  //
1021  // Get Psf candidate
1022  //
1023  std::shared_ptr<afwDetection::Psf::Image> kImage = psf.computeImage(afwGeom::PointD(x, y));
1024 
1025  //
1026  // Now find the proper sub-Image
1027  //
1028  afwGeom::BoxI bbox = kImage->getBBox();
1029 
1030  std::shared_ptr<MaskedImageT> subData(new MaskedImageT(*data, bbox, afwImage::PARENT, false)); // shallow copy
1031  //
1032  // Now we've got both; find the PSF's amplitude
1033  //
1034  double lambda = 0.0; // floor for variance is lambda*data
1035  try {
1036  double chi2; // chi^2 for fit
1037  double amp; // estimate of amplitude of model at this point
1038 
1039  if (std::isnan(psfFlux)) {
1040  std::pair<double, double> result = fitKernel(*kImage, *subData, lambda, true);
1041  chi2 = result.first; // chi^2 for fit
1042  amp = result.second; // estimate of amplitude of model at this point
1043  } else {
1044  chi2 = std::numeric_limits<double>::quiet_NaN();
1045  amp = psfFlux/afwMath::makeStatistics(*kImage, afwMath::SUM).getValue();
1046  }
1047  //
1048  // Convert kImage to the proper type so that I can subtract it.
1049  //
1050  std::shared_ptr<typename MaskedImageT::Image>
1051  kImageF(new typename MaskedImageT::Image(*kImage, true)); // of data's type
1052 
1053  *kImageF *= amp;
1054  *subData->getImage() -= *kImageF;
1055 
1056  return chi2;
1057  } catch(lsst::pex::exceptions::RangeError &e) {
1058  LSST_EXCEPT_ADD(e, (boost::format("Object at (%.2f, %.2f)") % x % y).str());
1059  throw e;
1060  }
1061 }
1062 
1063 /************************************************************************************************************/
1069 template<typename Image>
1070 std::pair<std::vector<double>, afwMath::KernelList>
1072  afwMath::LinearCombinationKernel const& kernel,
1073  Image const& image,
1074  afwGeom::Point2D const& pos
1075  )
1076 {
1078 
1079  afwMath::KernelList kernels = kernel.getKernelList(); // the Kernels that kernel adds together
1080  int const nKernel = kernels.size();
1081 
1082  if (nKernel == 0) {
1083  throw LSST_EXCEPT(lsst::pex::exceptions::LengthError,
1084  "Your kernel must have at least one component");
1085  }
1086 
1087  /*
1088  * Go through all the kernels, get a copy centered at the desired sub-pixel position, and then
1089  * extract a subImage from the parent image at the same place
1090  */
1091  std::vector<std::shared_ptr<KernelT>> kernelImages = offsetKernel<KernelT>(kernel, pos[0], pos[1]);
1092  afwGeom::BoxI bbox(kernelImages[0]->getBBox());
1093  Image const& subImage(Image(image, bbox, afwImage::PARENT, false)); // shallow copy
1094 
1095  /*
1096  * Solve the linear problem subImage = sum x_i K_i + epsilon; we solve this for x_i by constructing the
1097  * normal equations, A x = b
1098  */
1099  Eigen::MatrixXd A(nKernel, nKernel);
1100  Eigen::VectorXd b(nKernel);
1101 
1102  for (int i = 0; i != nKernel; ++i) {
1103  b(i) = afwImage::innerProduct(*kernelImages[i], *subImage.getImage());
1104 
1105  for (int j = i; j != nKernel; ++j) {
1106  A(i, j) = A(j, i) = afwImage::innerProduct(*kernelImages[i], *kernelImages[j]);
1107  }
1108  }
1109  Eigen::VectorXd x(nKernel);
1110 
1111  if (nKernel == 1) {
1112  x(0) = b(0)/A(0, 0);
1113  } else {
1114  x = A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
1115  }
1116 
1117  // the XY0() point of the shifted Kernel basis functions
1118  int const x0 = kernelImages[0]->getX0(), y0 = kernelImages[0]->getY0();
1119 
1120  afwMath::KernelList newKernels(nKernel);
1121  std::vector<double> params(nKernel);
1122  for (int i = 0; i != nKernel; ++i) {
1123  std::shared_ptr<afwMath::Kernel> newKernel(new afwMath::FixedKernel(*kernelImages[i]));
1124  newKernel->setCtrX(x0 + static_cast<int>(newKernel->getWidth()/2));
1125  newKernel->setCtrY(y0 + static_cast<int>(newKernel->getHeight()/2));
1126 
1127  params[i] = x[i];
1128  newKernels[i] = newKernel;
1129  }
1130 
1131  return std::make_pair(params, newKernels);
1132 }
1133 
1134 
1135 /************************************************************************************************************/
1141 template<typename Image>
1142 std::pair<std::shared_ptr<afwMath::Kernel>, std::pair<double, double> >
1144  afwMath::LinearCombinationKernel const& kernel,
1145  Image const& image,
1146  afwGeom::Point2D const& pos
1147  )
1148 {
1149  std::pair<std::vector<double>, afwMath::KernelList> const fit =
1150  fitKernelParamsToImage(kernel, image, pos);
1151  std::vector<double> params = fit.first;
1152  afwMath::KernelList kernels = fit.second;
1153  int const nKernel = params.size();
1154  assert(kernels.size() == static_cast<unsigned int>(nKernel));
1155 
1156  double amp = 0.0;
1157  for (int i = 0; i != nKernel; ++i) {
1158  std::shared_ptr<afwMath::Kernel> base = kernels[i];
1159  std::shared_ptr<afwMath::FixedKernel> k = std::static_pointer_cast<afwMath::FixedKernel>(base);
1160  amp += params[i] * k->getSum();
1161  }
1162 
1163  std::shared_ptr<afwMath::Kernel> outputKernel(new afwMath::LinearCombinationKernel(kernels, params));
1164  double chisq = 0.0;
1165  outputKernel->setCtrX(kernels[0]->getCtrX());
1166  outputKernel->setCtrY(kernels[0]->getCtrY());
1167 
1168  return std::make_pair(outputKernel, std::make_pair(amp, chisq));
1169 }
1170 
1171 
1172 /************************************************************************************************************/
1173 //
1174 // Explicit instantiations
1175 //
1177  typedef float Pixel;
1178 
1179  template
1180  std::pair<std::shared_ptr<afwMath::LinearCombinationKernel>, std::vector<double> >
1181  createKernelFromPsfCandidates<Pixel>(afwMath::SpatialCellSet const&, afwGeom::Extent2I const&,
1182  afwGeom::Point2I const&, int const, int const, int const,
1183  int const, bool const, int const);
1184  template
1185  int countPsfCandidates<Pixel>(afwMath::SpatialCellSet const&, int const);
1186 
1187  template
1188  std::pair<bool, double>
1189  fitSpatialKernelFromPsfCandidates<Pixel>(afwMath::Kernel *, afwMath::SpatialCellSet const&,
1190  int const, double const, double const);
1191  template
1192  std::pair<bool, double>
1193  fitSpatialKernelFromPsfCandidates<Pixel>(afwMath::Kernel *, afwMath::SpatialCellSet const&, bool const,
1194  int const, double const, double const);
1195 
1196  template
1197  double subtractPsf(afwDetection::Psf const&, afwImage::MaskedImage<float> *, double, double, double);
1198 
1199  template
1200  std::pair<std::vector<double>, afwMath::KernelList>
1201  fitKernelParamsToImage(afwMath::LinearCombinationKernel const&,
1202  afwImage::MaskedImage<Pixel> const&, afwGeom::Point2D const&);
1203 
1204  template
1205  std::pair<std::shared_ptr<afwMath::Kernel>, std::pair<double, double> >
1206  fitKernelToImage(afwMath::LinearCombinationKernel const&,
1207  afwImage::MaskedImage<Pixel> const&, afwGeom::Point2D const&);
1209 }}}
Class used by SpatialCell for spatial PSF fittig.
boost::shared_ptr< afw::image::MaskedImage< PixelT > > getOffsetImage(std::string const algorithm, unsigned int buffer) const
Return an offset version of the image of the source.
double operator()(const std::vector< double > &coeffs) const
void setAmplitude(double amplitude)
Set the best-fit amplitude.
Definition: PsfCandidate.h:116
void processCandidate(afwMath::SpatialCellCandidate *candidate)
int countPsfCandidates(lsst::afw::math::SpatialCellSet const &psfCells, int const nStarPerCell=-1)
boost::shared_ptr< afw::image::MaskedImage< PixelT > const > getMaskedImage() const
Return the image at the position of the Source, without any sub-pixel shifts to put the centre of the...
double subtractPsf(lsst::afw::detection::Psf const &psf, ImageT *data, double x, double y, double psfFlux=std::numeric_limits< double >::quiet_NaN())
MinimizeChi2(evalChi2Visitor< PixelT > &chi2Visitor, afwMath::Kernel *kernel, afwMath::SpatialCellSet const &psfCells, int nStarPerCell, int nComponents, int nSpatialParams)
STL namespace.
virtual void analyze()
Generate eigenimages that are normalised and background-subtracted.
Definition: ImagePca.cc:41
evalChi2Visitor(afwMath::Kernel const &kernel, double lambda)
std::pair< std::vector< double >, lsst::afw::math::KernelList > fitKernelParamsToImage(lsst::afw::math::LinearCombinationKernel const &kernel, Image const &image, lsst::afw::geom::Point2D const &pos)
std::pair< bool, double > fitSpatialKernelFromPsfCandidates(lsst::afw::math::Kernel *kernel, lsst::afw::math::SpatialCellSet const &psfCells, int const nStarPerCell=-1, double const tolerance=1e-5, double const lambda=0.0)
boost::shared_ptr< afw::table::SourceRecord > getSource() const
Return the original Source.
Definition: PsfCandidate.h:110
double getVar() const
Return the variance in use when fitting this object.
Definition: PsfCandidate.h:119
double getAmplitude() const
Return the best-fit amplitude.
Definition: PsfCandidate.h:113
A class to pass around to all our PsfCandidates to evaluate the PSF fit&#39;s X^2.
std::pair< std::shared_ptr< lsst::afw::math::LinearCombinationKernel >, std::vector< double > > createKernelFromPsfCandidates(lsst::afw::math::SpatialCellSet const &psfCells, lsst::afw::geom::Extent2I const &dims, lsst::afw::geom::Point2I const &xy0, int const nEigenComponents, int const spatialOrder, int const ksize, int const nStarPerCell=-1, bool const constantWeight=true, int const border=3)
Class used by SpatialCell for spatial PSF fittig.
Class for doing PCA on PSF stars.
void setSpatialParameters(afwMath::Kernel *kernel, std::vector< double > const &coeffs)
Fit a Kernel&#39;s spatial variability from a set of stars.
double Up() const
Error definition of the function.
afw::table::Key< double > b
std::pair< std::shared_ptr< lsst::afw::math::Kernel >, std::pair< double, double > > fitKernelToImage(lsst::afw::math::LinearCombinationKernel const &kernel, Image const &image, lsst::afw::geom::Point2D const &pos)
Class stored in SpatialCells for spatial Psf fitting.
Definition: PsfCandidate.h:56