LSSTApplications  10.0+286,10.0+36,10.0+46,10.0-2-g4f67435,10.1+152,10.1+37,11.0,11.0+1,11.0-1-g47edd16,11.0-1-g60db491,11.0-1-g7418c06,11.0-2-g04d2804,11.0-2-g68503cd,11.0-2-g818369d,11.0-2-gb8b8ce7
LSSTDataManagementBasePackage
KernelSolution.cc
Go to the documentation of this file.
1 // -*- lsst-c++ -*-
11 #include <iterator>
12 #include <cmath>
13 #include <algorithm>
14 #include <limits>
15 
16 #include "boost/shared_ptr.hpp"
17 #include "boost/timer.hpp"
18 
19 #include "Eigen/Core"
20 #include "Eigen/Cholesky"
21 #include "Eigen/QR"
22 #include "Eigen/LU"
23 #include "Eigen/Eigenvalues"
24 #include "Eigen/SVD"
25 
26 #include "lsst/afw/math.h"
27 #include "lsst/afw/geom.h"
28 #include "lsst/afw/image.h"
29 #include "lsst/afw/detection.h"
30 #include "lsst/afw/detection/FootprintArray.cc"
32 #include "lsst/pex/logging/Trace.h"
33 
36 
37 #include "ndarray.h"
38 #include "ndarray/eigen.h"
39 
40 #define DEBUG_MATRIX 0
41 #define DEBUG_MATRIX2 0
42 
43 namespace afwDet = lsst::afw::detection;
44 namespace afwMath = lsst::afw::math;
45 namespace afwGeom = lsst::afw::geom;
46 namespace afwImage = lsst::afw::image;
47 namespace pexLog = lsst::pex::logging;
48 namespace pexExcept = lsst::pex::exceptions;
49 
50 namespace lsst {
51 namespace ip {
52 namespace diffim {
53 
54  /* Unique identifier for solution */
56 
58  boost::shared_ptr<Eigen::MatrixXd> mMat,
59  boost::shared_ptr<Eigen::VectorXd> bVec,
60  bool fitForBackground
61  ) :
62  _id(++_SolutionId),
63  _mMat(mMat),
64  _bVec(bVec),
65  _aVec(),
66  _solvedBy(NONE),
67  _fitForBackground(fitForBackground)
68  {};
69 
71  bool fitForBackground
72  ) :
73  _id(++_SolutionId),
74  _mMat(),
75  _bVec(),
76  _aVec(),
77  _solvedBy(NONE),
78  _fitForBackground(fitForBackground)
79  {};
80 
82  _id(++_SolutionId),
83  _mMat(),
84  _bVec(),
85  _aVec(),
86  _solvedBy(NONE),
87  _fitForBackground(true)
88  {};
89 
91  solve(*_mMat, *_bVec);
92  }
93 
95  return getConditionNumber(*_mMat, conditionType);
96  }
97 
98  double KernelSolution::getConditionNumber(Eigen::MatrixXd mMat,
99  ConditionNumberType conditionType) {
100  switch (conditionType) {
101  case EIGENVALUE:
102  {
103  Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(mMat);
104  Eigen::VectorXd eValues = eVecValues.eigenvalues();
105  double eMax = eValues.maxCoeff();
106  double eMin = eValues.minCoeff();
107  pexLog::TTrace<5>("lsst.ip.diffim.KernelSolution.getConditionNumber",
108  "EIGENVALUE eMax / eMin = %.3e", eMax / eMin);
109  return (eMax / eMin);
110  break;
111  }
112  case SVD:
113  {
114  Eigen::VectorXd sValues = mMat.jacobiSvd().singularValues();
115  double sMax = sValues.maxCoeff();
116  double sMin = sValues.minCoeff();
117  pexLog::TTrace<5>("lsst.ip.diffim.KernelSolution.getConditionNumber",
118  "SVD eMax / eMin = %.3e", sMax / sMin);
119  return (sMax / sMin);
120  break;
121  }
122  default:
123  {
124  throw LSST_EXCEPT(pexExcept::InvalidParameterError,
125  "Undefined ConditionNumberType : only EIGENVALUE, SVD allowed.");
126  break;
127  }
128  }
129  }
130 
131  void KernelSolution::solve(Eigen::MatrixXd mMat,
132  Eigen::VectorXd bVec) {
133 
134  if (DEBUG_MATRIX) {
135  std::cout << "M " << std::endl;
136  std::cout << mMat << std::endl;
137  std::cout << "B " << std::endl;
138  std::cout << bVec << std::endl;
139  }
140 
141  Eigen::VectorXd aVec = Eigen::VectorXd::Zero(bVec.size());
142 
143  boost::timer t;
144  t.restart();
145 
146  pexLog::TTrace<4>("lsst.ip.difim.KernelSolution.solve",
147  "Solving for kernel");
148  _solvedBy = LU;
149  Eigen::FullPivLU<Eigen::MatrixXd> lu(mMat);
150  if (lu.isInvertible()) {
151  aVec = lu.solve(bVec);
152  } else {
153  pexLog::TTrace<5>("lsst.ip.diffim.KernelSolution.solve",
154  "Unable to determine kernel via LU");
155  /* LAST RESORT */
156  try {
157 
158  _solvedBy = EIGENVECTOR;
159  Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(mMat);
160  Eigen::MatrixXd const& rMat = eVecValues.eigenvectors();
161  Eigen::VectorXd eValues = eVecValues.eigenvalues();
162 
163  for (int i = 0; i != eValues.rows(); ++i) {
164  if (eValues(i) != 0.0) {
165  eValues(i) = 1.0/eValues(i);
166  }
167  }
168 
169  aVec = rMat * eValues.asDiagonal() * rMat.transpose() * bVec;
170  } catch (pexExcept::Exception& e) {
171 
172  _solvedBy = NONE;
173  pexLog::TTrace<5>("lsst.ip.diffim.KernelSolution.solve",
174  "Unable to determine kernel via eigen-values");
175 
176  throw LSST_EXCEPT(pexExcept::Exception, "Unable to determine kernel solution");
177  }
178  }
179 
180  double time = t.elapsed();
181  pexLog::TTrace<5>("lsst.ip.diffim.KernelSolution.solve",
182  "Compute time for matrix math : %.2f s", time);
183 
184  if (DEBUG_MATRIX) {
185  std::cout << "A " << std::endl;
186  std::cout << aVec << std::endl;
187  }
188 
189  _aVec = boost::shared_ptr<Eigen::VectorXd>(new Eigen::VectorXd(aVec));
190  }
191 
192  /*******************************************************************************************************/
193 
194  template <typename InputT>
196  lsst::afw::math::KernelList const& basisList,
197  bool fitForBackground
198  )
199  :
200  KernelSolution(fitForBackground),
201  _cMat(),
202  _iVec(),
203  _ivVec(),
204  _kernel(),
205  _background(0.0),
206  _kSum(0.0)
207  {
208  std::vector<double> kValues(basisList.size());
209  _kernel = boost::shared_ptr<afwMath::Kernel>(
210  new afwMath::LinearCombinationKernel(basisList, kValues)
211  );
212  };
213 
214  template <typename InputT>
216  if (_solvedBy == KernelSolution::NONE) {
217  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return solution");
218  }
219  return _kernel;
220  }
221 
222  template <typename InputT>
224  if (_solvedBy == KernelSolution::NONE) {
225  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return image");
226  }
229  );
230  (void)_kernel->computeImage(*image, false);
231  return image;
232  }
233 
234  template <typename InputT>
236  if (_solvedBy == KernelSolution::NONE) {
237  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return background");
238  }
239  return _background;
240  }
241 
242  template <typename InputT>
244  if (_solvedBy == KernelSolution::NONE) {
245  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return ksum");
246  }
247  return _kSum;
248  }
249 
250  template <typename InputT>
251  std::pair<boost::shared_ptr<lsst::afw::math::Kernel>, double>
253  if (_solvedBy == KernelSolution::NONE) {
254  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return solution");
255  }
256 
257  return std::make_pair(_kernel, _background);
258  }
259 
260  template <typename InputT>
262  lsst::afw::image::Image<InputT> const &templateImage,
263  lsst::afw::image::Image<InputT> const &scienceImage,
265  ) {
266 
267  afwMath::Statistics varStats = afwMath::makeStatistics(varianceEstimate, afwMath::MIN);
268  if (varStats.getValue(afwMath::MIN) < 0.0) {
270  "Error: variance less than 0.0");
271  }
272  if (varStats.getValue(afwMath::MIN) == 0.0) {
274  "Error: variance equals 0.0, cannot inverse variance weight");
275  }
276 
277  lsst::afw::math::KernelList basisList =
278  boost::dynamic_pointer_cast<afwMath::LinearCombinationKernel>(_kernel)->getKernelList();
279 
280  unsigned int const nKernelParameters = basisList.size();
281  unsigned int const nBackgroundParameters = _fitForBackground ? 1 : 0;
282  unsigned int const nParameters = nKernelParameters + nBackgroundParameters;
283 
284  std::vector<boost::shared_ptr<afwMath::Kernel> >::const_iterator kiter = basisList.begin();
285 
286  /* Ignore buffers around edge of convolved images :
287  *
288  * If the kernel has width 5, it has center pixel 2. The first good pixel
289  * is the (5-2)=3rd pixel, which is array index 2, and ends up being the
290  * index of the central pixel.
291  *
292  * You also have a buffer of unusable pixels on the other side, numbered
293  * width-center-1. The last good usable pixel is N-width+center+1.
294  *
295  * Example : the kernel is width = 5, center = 2
296  *
297  * ---|---|-c-|---|---|
298  *
299  * the image is width = N
300  * convolve this with the kernel, and you get
301  *
302  * |-x-|-x-|-g-|---|---| ... |---|---|-g-|-x-|-x-|
303  *
304  * g = first/last good pixel
305  * x = bad
306  *
307  * the first good pixel is the array index that has the value "center", 2
308  * the last good pixel has array index N-(5-2)+1
309  * eg. if N = 100, you want to use up to index 97
310  * 100-3+1 = 98, and the loops use i < 98, meaning the last
311  * index you address is 97.
312  */
313 
314  /* NOTE - we are accessing particular elements of Eigen arrays using
315  these coordinates, therefore they need to be in LOCAL coordinates.
316  This was written before ndarray unification.
317  */
318  afwGeom::Box2I goodBBox = (*kiter)->shrinkBBox(templateImage.getBBox(afwImage::LOCAL));
319  unsigned int const startCol = goodBBox.getMinX();
320  unsigned int const startRow = goodBBox.getMinY();
321  // endCol/Row is one past the index of the last good col/row
322  unsigned int endCol = goodBBox.getMaxX() + 1;
323  unsigned int endRow = goodBBox.getMaxY() + 1;
324 
325  boost::timer t;
326  t.restart();
327 
328  /* Eigen representation of input images; only the pixels that are unconvolved in cimage below */
329  Eigen::MatrixXd eigenTemplate = imageToEigenMatrix(templateImage).block(startRow,
330  startCol,
331  endRow-startRow,
332  endCol-startCol);
333  Eigen::MatrixXd eigenScience = imageToEigenMatrix(scienceImage).block(startRow,
334  startCol,
335  endRow-startRow,
336  endCol-startCol);
337  Eigen::MatrixXd eigeniVariance = imageToEigenMatrix(varianceEstimate).block(
338  startRow, startCol, endRow-startRow, endCol-startCol
339  ).array().inverse().matrix();
340 
341  /* Resize into 1-D for later usage */
342  eigenTemplate.resize(eigenTemplate.rows()*eigenTemplate.cols(), 1);
343  eigenScience.resize(eigenScience.rows()*eigenScience.cols(), 1);
344  eigeniVariance.resize(eigeniVariance.rows()*eigeniVariance.cols(), 1);
345 
346  /* Holds image convolved with basis function */
347  afwImage::Image<PixelT> cimage(templateImage.getDimensions());
348 
349  /* Holds eigen representation of image convolved with all basis functions */
350  std::vector<boost::shared_ptr<Eigen::MatrixXd> > convolvedEigenList(nKernelParameters);
351 
352  /* Iterators over convolved image list and basis list */
353  typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiter =
354  convolvedEigenList.begin();
355  /* Create C_i in the formalism of Alard & Lupton */
356  for (kiter = basisList.begin(); kiter != basisList.end(); ++kiter, ++eiter) {
357  afwMath::convolve(cimage, templateImage, **kiter, false); /* cimage stores convolved image */
358 
359  boost::shared_ptr<Eigen::MatrixXd> cMat (
360  new Eigen::MatrixXd(imageToEigenMatrix(cimage).block(startRow,
361  startCol,
362  endRow-startRow,
363  endCol-startCol))
364  );
365  cMat->resize(cMat->rows()*cMat->cols(), 1);
366  *eiter = cMat;
367 
368  }
369 
370  double time = t.elapsed();
371  pexLog::TTrace<5>("lsst.ip.diffim.StaticKernelSolution.build",
372  "Total compute time to do basis convolutions : %.2f s", time);
373  t.restart();
374 
375  /*
376  Load matrix with all values from convolvedEigenList : all images
377  (eigeniVariance, convolvedEigenList) must be the same size
378  */
379  Eigen::MatrixXd cMat(eigenTemplate.col(0).size(), nParameters);
380  typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiterj =
381  convolvedEigenList.begin();
382  typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiterE =
383  convolvedEigenList.end();
384  for (unsigned int kidxj = 0; eiterj != eiterE; eiterj++, kidxj++) {
385  cMat.col(kidxj) = (*eiterj)->col(0);
386  }
387  /* Treat the last "image" as all 1's to do the background calculation. */
388  if (_fitForBackground)
389  cMat.col(nParameters-1).fill(1.);
390 
391  _cMat.reset(new Eigen::MatrixXd(cMat));
392  _ivVec.reset(new Eigen::VectorXd(eigeniVariance.col(0)));
393  _iVec.reset(new Eigen::VectorXd(eigenScience.col(0)));
394 
395  /* Make these outside of solve() so I can check condition number */
396  _mMat.reset(new Eigen::MatrixXd((*_cMat).transpose() * ((*_ivVec).asDiagonal() * (*_cMat))));
397  _bVec.reset(new Eigen::VectorXd((*_cMat).transpose() * ((*_ivVec).asDiagonal() * (*_iVec))));
398  }
399 
400  template <typename InputT>
402  pexLog::TTrace<5>("lsst.ip.diffim.StaticKernelSolution.solve",
403  "mMat is %d x %d; bVec is %d; cMat is %d x %d; vVec is %d; iVec is %d",
404  (*_mMat).rows(), (*_mMat).cols(), (*_bVec).size(),
405  (*_cMat).rows(), (*_cMat).cols(), (*_ivVec).size(), (*_iVec).size());
406 
407  /* If I put this here I can't check for condition number before solving */
408  /*
409  _mMat.reset(new Eigen::MatrixXd((*_cMat).transpose() * ((*_ivVec).asDiagonal() * (*_cMat))));
410  _bVec.reset(new Eigen::VectorXd((*_cMat).transpose() * ((*_ivVec).asDiagonal() * (*_iVec))));
411  */
412 
413  if (DEBUG_MATRIX) {
414  std::cout << "C" << std::endl;
415  std::cout << (*_cMat) << std::endl;
416  std::cout << "iV" << std::endl;
417  std::cout << (*_ivVec) << std::endl;
418  std::cout << "I" << std::endl;
419  std::cout << (*_iVec) << std::endl;
420  }
421 
422  try {
424  } catch (pexExcept::Exception &e) {
425  LSST_EXCEPT_ADD(e, "Unable to solve static kernel matrix");
426  throw e;
427  }
428  /* Turn matrices into _kernel and _background */
429  _setKernel();
430  }
431 
432  template <typename InputT>
434  if (_solvedBy == KernelSolution::NONE) {
435  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot make solution");
436  }
437 
438  unsigned int const nParameters = _aVec->size();
439  unsigned int const nBackgroundParameters = _fitForBackground ? 1 : 0;
440  unsigned int const nKernelParameters =
441  boost::dynamic_pointer_cast<afwMath::LinearCombinationKernel>(_kernel)->getKernelList().size();
442  if (nParameters != (nKernelParameters + nBackgroundParameters))
443  throw LSST_EXCEPT(pexExcept::Exception, "Mismatched sizes in kernel solution");
444 
445  /* Fill in the kernel results */
446  std::vector<double> kValues(nKernelParameters);
447  for (unsigned int idx = 0; idx < nKernelParameters; idx++) {
448  if (std::isnan((*_aVec)(idx))) {
450  str(boost::format("Unable to determine kernel solution %d (nan)") % idx));
451  }
452  kValues[idx] = (*_aVec)(idx);
453  }
454  _kernel->setKernelParameters(kValues);
455 
458  );
459  _kSum = _kernel->computeImage(*image, false);
460 
461  if (_fitForBackground) {
462  if (std::isnan((*_aVec)(nParameters-1))) {
464  str(boost::format("Unable to determine background solution %d (nan)") %
465  (nParameters-1)));
466  }
467  _background = (*_aVec)(nParameters-1);
468  }
469  }
470 
471 
472  template <typename InputT>
474  throw LSST_EXCEPT(pexExcept::Exception, "Uncertainty calculation not supported");
475 
476  /* Estimate of parameter uncertainties comes from the inverse of the
477  * covariance matrix (noise spectrum).
478  * N.R. 15.4.8 to 15.4.15
479  *
480  * Since this is a linear problem no need to use Fisher matrix
481  * N.R. 15.5.8
482  *
483  * Although I might be able to take advantage of the solution above.
484  * Since this now works and is not the rate limiting step, keep as-is for DC3a.
485  *
486  * Use Cholesky decomposition again.
487  * Cholkesy:
488  * Cov = L L^t
489  * Cov^(-1) = (L L^t)^(-1)
490  * = (L^T)^-1 L^(-1)
491  *
492  * Code would be:
493  *
494  * Eigen::MatrixXd cov = (*_mMat).transpose() * (*_mMat);
495  * Eigen::LLT<Eigen::MatrixXd> llt = cov.llt();
496  * Eigen::MatrixXd error2 = llt.matrixL().transpose().inverse()*llt.matrixL().inverse();
497  */
498  }
499 
500  /*******************************************************************************************************/
501 
502  template <typename InputT>
504  lsst::afw::math::KernelList const& basisList,
505  bool fitForBackground
506  ) :
507  StaticKernelSolution<InputT>(basisList, fitForBackground)
508  {};
509 
510  template <typename InputT>
512  lsst::afw::image::Image<InputT> const &templateImage,
513  lsst::afw::image::Image<InputT> const &scienceImage,
516  ) {
517 
518  afwMath::Statistics varStats = afwMath::makeStatistics(varianceEstimate, afwMath::MIN);
519  if (varStats.getValue(afwMath::MIN) < 0.0) {
521  "Error: variance less than 0.0");
522  }
523  if (varStats.getValue(afwMath::MIN) == 0.0) {
525  "Error: variance equals 0.0, cannot inverse variance weight");
526  }
527 
528  /* Full footprint of all input images */
529  afwDet::Footprint::Ptr fullFp(new afwDet::Footprint(templateImage.getBBox()));
530 
531  afwMath::KernelList basisList =
532  boost::dynamic_pointer_cast<afwMath::LinearCombinationKernel>(this->_kernel)->getKernelList();
533  std::vector<boost::shared_ptr<afwMath::Kernel> >::const_iterator kiter = basisList.begin();
534 
535  /* Only BAD pixels marked in this mask */
536  afwImage::MaskPixel bitMask =
541 
542  /* Create a Footprint that contains all the masked pixels set above */
544  afwDet::FootprintSet maskFpSet(pixelMask, threshold, true);
545 
546  /* And spread it by the kernel half width */
547  int growPix = (*kiter)->getCtr().getX();
548  afwDet::FootprintSet maskedFpSetGrown(maskFpSet, growPix, true);
549 
550 #if 0
551  for (typename afwDet::FootprintSet::FootprintList::iterator
552  ptr = maskedFpSetGrown.getFootprints()->begin(),
553  end = maskedFpSetGrown.getFootprints()->end();
554  ptr != end;
555  ++ptr) {
556 
557  afwDet::setMaskFromFootprint(finalMask,
558  (**ptr),
560  }
561 #endif
562 
563  afwImage::Mask<afwImage::MaskPixel> finalMask(pixelMask.getDimensions());
565  *(maskedFpSetGrown.getFootprints()),
567  pixelMask.writeFits("pixelmask.fits");
568  finalMask.writeFits("finalmask.fits");
569 
570 
571  ndarray::Array<int, 1, 1> maskArray =
572  ndarray::allocate(ndarray::makeVector(fullFp->getArea()));
573  afwDet::flattenArray(*fullFp, finalMask.getArray(),
574  maskArray, templateImage.getXY0()); /* Need to fake the XY0 */
575  ndarray::EigenView<int, 1, 1> maskEigen(maskArray);
576 
577  ndarray::Array<InputT, 1, 1> arrayTemplate =
578  ndarray::allocate(ndarray::makeVector(fullFp->getArea()));
579  afwDet::flattenArray(*fullFp, templateImage.getArray(),
580  arrayTemplate, templateImage.getXY0());
581  ndarray::EigenView<InputT, 1, 1> eigenTemplate0(arrayTemplate);
582 
583  ndarray::Array<InputT, 1, 1> arrayScience =
584  ndarray::allocate(ndarray::makeVector(fullFp->getArea()));
585  afwDet::flattenArray(*fullFp, scienceImage.getArray(),
586  arrayScience, scienceImage.getXY0());
587  ndarray::EigenView<InputT, 1, 1> eigenScience0(arrayScience);
588 
590  ndarray::allocate(ndarray::makeVector(fullFp->getArea()));
591  afwDet::flattenArray(*fullFp, varianceEstimate.getArray(),
592  arrayVariance, varianceEstimate.getXY0());
593  ndarray::EigenView<afwImage::VariancePixel, 1, 1> eigenVariance0(arrayVariance);
594 
595  int nGood = 0;
596  for (int i = 0; i < maskEigen.size(); i++) {
597  if (maskEigen(i) == 0.0)
598  nGood += 1;
599  }
600 
601  Eigen::VectorXd eigenTemplate(nGood);
602  Eigen::VectorXd eigenScience(nGood);
603  Eigen::VectorXd eigenVariance(nGood);
604  int nUsed = 0;
605  for (int i = 0; i < maskEigen.size(); i++) {
606  if (maskEigen(i) == 0.0) {
607  eigenTemplate(nUsed) = eigenTemplate0(i);
608  eigenScience(nUsed) = eigenScience0(i);
609  eigenVariance(nUsed) = eigenVariance0(i);
610  nUsed += 1;
611  }
612  }
613 
614 
615  boost::timer t;
616  t.restart();
617 
618  unsigned int const nKernelParameters = basisList.size();
619  unsigned int const nBackgroundParameters = this->_fitForBackground ? 1 : 0;
620  unsigned int const nParameters = nKernelParameters + nBackgroundParameters;
621 
622  /* Holds image convolved with basis function */
623  afwImage::Image<InputT> cimage(templateImage.getDimensions());
624 
625  /* Holds eigen representation of image convolved with all basis functions */
626  std::vector<boost::shared_ptr<Eigen::VectorXd> >
627  convolvedEigenList(nKernelParameters);
628 
629  /* Iterators over convolved image list and basis list */
630  typename std::vector<boost::shared_ptr<Eigen::VectorXd> >::iterator eiter =
631  convolvedEigenList.begin();
632 
633  /* Create C_i in the formalism of Alard & Lupton */
634  for (kiter = basisList.begin(); kiter != basisList.end(); ++kiter, ++eiter) {
635  afwMath::convolve(cimage, templateImage, **kiter, false); /* cimage stores convolved image */
636 
638  ndarray::allocate(ndarray::makeVector(fullFp->getArea()));
639  afwDet::flattenArray(*fullFp, cimage.getArray(),
640  arrayC, cimage.getXY0());
641  ndarray::EigenView<InputT, 1, 1> eigenC0(arrayC);
642 
643  Eigen::VectorXd eigenC(nGood);
644  int nUsed = 0;
645  for (int i = 0; i < maskEigen.size(); i++) {
646  if (maskEigen(i) == 0.0) {
647  eigenC(nUsed) = eigenC0(i);
648  nUsed += 1;
649  }
650  }
651 
652  boost::shared_ptr<Eigen::VectorXd>
653  eigenCptr (new Eigen::VectorXd(eigenC));
654 
655  *eiter = eigenCptr;
656  }
657  double time = t.elapsed();
658  pexLog::TTrace<5>("lsst.ip.diffim.StaticKernelSolution.buildWithMask",
659  "Total compute time to do basis convolutions : %.2f s", time);
660  t.restart();
661 
662  /* Load matrix with all convolved images */
663  Eigen::MatrixXd cMat(eigenTemplate.size(), nParameters);
664  typename std::vector<boost::shared_ptr<Eigen::VectorXd> >::iterator eiterj =
665  convolvedEigenList.begin();
666  typename std::vector<boost::shared_ptr<Eigen::VectorXd> >::iterator eiterE =
667  convolvedEigenList.end();
668  for (unsigned int kidxj = 0; eiterj != eiterE; eiterj++, kidxj++) {
669  cMat.block(0, kidxj, eigenTemplate.size(), 1) =
670  Eigen::MatrixXd(**eiterj).block(0, 0, eigenTemplate.size(), 1);
671  }
672  /* Treat the last "image" as all 1's to do the background calculation. */
673  if (this->_fitForBackground)
674  cMat.col(nParameters-1).fill(1.);
675 
676  this->_cMat.reset(new Eigen::MatrixXd(cMat));
677  //this->_ivVec.reset(new Eigen::VectorXd((eigenVariance.template cast<double>()).cwise().inverse()));
678  //this->_iVec.reset(new Eigen::VectorXd(eigenScience.template cast<double>()));
679  this->_ivVec.reset(new Eigen::VectorXd(eigenVariance.array().inverse().matrix()));
680  this->_iVec.reset(new Eigen::VectorXd(eigenScience));
681 
682  /* Make these outside of solve() so I can check condition number */
683  this->_mMat.reset(
684  new Eigen::MatrixXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_cMat))
685  );
686  this->_bVec.reset(
687  new Eigen::VectorXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_iVec))
688  );
689  }
690 
691 
692  template <typename InputT>
694  lsst::afw::image::Image<InputT> const &templateImage,
695  lsst::afw::image::Image<InputT> const &scienceImage,
698  ) {
699 
700  afwMath::Statistics varStats = afwMath::makeStatistics(varianceEstimate, afwMath::MIN);
701  if (varStats.getValue(afwMath::MIN) < 0.0) {
703  "Error: variance less than 0.0");
704  }
705  if (varStats.getValue(afwMath::MIN) == 0.0) {
707  "Error: variance equals 0.0, cannot inverse variance weight");
708  }
709 
710  lsst::afw::math::KernelList basisList =
711  boost::dynamic_pointer_cast<afwMath::LinearCombinationKernel>(this->_kernel)->getKernelList();
712  std::vector<boost::shared_ptr<afwMath::Kernel> >::const_iterator kiter = basisList.begin();
713 
714  /* Only BAD pixels marked in this mask */
716  afwImage::MaskPixel bitMask =
720  sMask &= bitMask;
721  /* TBD: Need to figure out a way to spread this; currently done in Python */
722 
723  unsigned int const nKernelParameters = basisList.size();
724  unsigned int const nBackgroundParameters = this->_fitForBackground ? 1 : 0;
725  unsigned int const nParameters = nKernelParameters + nBackgroundParameters;
726 
727  /* NOTE - we are accessing particular elements of Eigen arrays using
728  these coordinates, therefore they need to be in LOCAL coordinates.
729  This was written before ndarray unification.
730  */
731  /* Ignore known EDGE pixels for speed */
732  afwGeom::Box2I shrunkLocalBBox = (*kiter)->shrinkBBox(templateImage.getBBox(afwImage::LOCAL));
733  pexLog::TTrace<5>("lsst.ip.diffim.MaskedKernelSolution.build",
734  "Limits of good pixels after convolution: %d,%d -> %d,%d (local)",
735  shrunkLocalBBox.getMinX(), shrunkLocalBBox.getMinY(),
736  shrunkLocalBBox.getMaxX(), shrunkLocalBBox.getMaxY());
737 
738  /* Subimages are addressed in raw pixel coordinates */
739  unsigned int startCol = shrunkLocalBBox.getMinX();
740  unsigned int startRow = shrunkLocalBBox.getMinY();
741  unsigned int endCol = shrunkLocalBBox.getMaxX();
742  unsigned int endRow = shrunkLocalBBox.getMaxY();
743  /* Needed for Eigen block slicing operation */
744  endCol += 1;
745  endRow += 1;
746 
747  boost::timer t;
748  t.restart();
749 
750  /* Eigen representation of input images; only the pixels that are unconvolved in cimage below */
751  Eigen::MatrixXi eMask = maskToEigenMatrix(sMask).block(startRow,
752  startCol,
753  endRow-startRow,
754  endCol-startCol);
755 
756  Eigen::MatrixXd eigenTemplate = imageToEigenMatrix(templateImage).block(startRow,
757  startCol,
758  endRow-startRow,
759  endCol-startCol);
760  Eigen::MatrixXd eigenScience = imageToEigenMatrix(scienceImage).block(startRow,
761  startCol,
762  endRow-startRow,
763  endCol-startCol);
764  Eigen::MatrixXd eigeniVariance = imageToEigenMatrix(varianceEstimate).block(
765  startRow, startCol, endRow-startRow, endCol-startCol
766  ).array().inverse().matrix();
767 
768  /* Resize into 1-D for later usage */
769  eMask.resize(eMask.rows()*eMask.cols(), 1);
770  eigenTemplate.resize(eigenTemplate.rows()*eigenTemplate.cols(), 1);
771  eigenScience.resize(eigenScience.rows()*eigenScience.cols(), 1);
772  eigeniVariance.resize(eigeniVariance.rows()*eigeniVariance.cols(), 1);
773 
774  /* Do the masking, slowly for now... */
775  Eigen::MatrixXd maskedEigenTemplate(eigenTemplate.rows(), 1);
776  Eigen::MatrixXd maskedEigenScience(eigenScience.rows(), 1);
777  Eigen::MatrixXd maskedEigeniVariance(eigeniVariance.rows(), 1);
778  int nGood = 0;
779  for (int i = 0; i < eMask.rows(); i++) {
780  if (eMask(i, 0) == 0) {
781  maskedEigenTemplate(nGood, 0) = eigenTemplate(i, 0);
782  maskedEigenScience(nGood, 0) = eigenScience(i, 0);
783  maskedEigeniVariance(nGood, 0) = eigeniVariance(i, 0);
784  nGood += 1;
785  }
786  }
787  /* Can't resize() since all values are lost; use blocks */
788  eigenTemplate = maskedEigenTemplate.block(0, 0, nGood, 1);
789  eigenScience = maskedEigenScience.block(0, 0, nGood, 1);
790  eigeniVariance = maskedEigeniVariance.block(0, 0, nGood, 1);
791 
792 
793  /* Holds image convolved with basis function */
794  afwImage::Image<InputT> cimage(templateImage.getDimensions());
795 
796  /* Holds eigen representation of image convolved with all basis functions */
797  std::vector<boost::shared_ptr<Eigen::MatrixXd> > convolvedEigenList(nKernelParameters);
798 
799  /* Iterators over convolved image list and basis list */
800  typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiter =
801  convolvedEigenList.begin();
802  /* Create C_i in the formalism of Alard & Lupton */
803  for (kiter = basisList.begin(); kiter != basisList.end(); ++kiter, ++eiter) {
804  afwMath::convolve(cimage, templateImage, **kiter, false); /* cimage stores convolved image */
805 
806  Eigen::MatrixXd cMat = imageToEigenMatrix(cimage).block(startRow,
807  startCol,
808  endRow-startRow,
809  endCol-startCol);
810  cMat.resize(cMat.rows() * cMat.cols(), 1);
811 
812  /* Do masking */
813  Eigen::MatrixXd maskedcMat(cMat.rows(), 1);
814  int nGood = 0;
815  for (int i = 0; i < eMask.rows(); i++) {
816  if (eMask(i, 0) == 0) {
817  maskedcMat(nGood, 0) = cMat(i, 0);
818  nGood += 1;
819  }
820  }
821  cMat = maskedcMat.block(0, 0, nGood, 1);
822  boost::shared_ptr<Eigen::MatrixXd> cMatPtr (new Eigen::MatrixXd(cMat));
823  *eiter = cMatPtr;
824  }
825 
826  double time = t.elapsed();
827  pexLog::TTrace<5>("lsst.ip.diffim.StaticKernelSolution.build",
828  "Total compute time to do basis convolutions : %.2f s", time);
829  t.restart();
830 
831  /*
832  Load matrix with all values from convolvedEigenList : all images
833  (eigeniVariance, convolvedEigenList) must be the same size
834  */
835  Eigen::MatrixXd cMat(eigenTemplate.col(0).size(), nParameters);
836  typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiterj =
837  convolvedEigenList.begin();
838  typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiterE =
839  convolvedEigenList.end();
840  for (unsigned int kidxj = 0; eiterj != eiterE; eiterj++, kidxj++) {
841  cMat.col(kidxj) = (*eiterj)->col(0);
842  }
843  /* Treat the last "image" as all 1's to do the background calculation. */
844  if (this->_fitForBackground)
845  cMat.col(nParameters-1).fill(1.);
846 
847  this->_cMat.reset(new Eigen::MatrixXd(cMat));
848  this->_ivVec.reset(new Eigen::VectorXd(eigeniVariance.col(0)));
849  this->_iVec.reset(new Eigen::VectorXd(eigenScience.col(0)));
850 
851  /* Make these outside of solve() so I can check condition number */
852  this->_mMat.reset(
853  new Eigen::MatrixXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_cMat))
854  );
855  this->_bVec.reset(
856  new Eigen::VectorXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_iVec))
857  );
858 
859  }
860 
861  /* NOTE - this was written before the ndarray unification. I am rewriting
862  buildSingleMask to use the ndarray stuff */
863  template <typename InputT>
865  lsst::afw::image::Image<InputT> const &templateImage,
866  lsst::afw::image::Image<InputT> const &scienceImage,
868  lsst::afw::geom::Box2I maskBox
869  ) {
870 
871  afwMath::Statistics varStats = afwMath::makeStatistics(varianceEstimate, afwMath::MIN);
872  if (varStats.getValue(afwMath::MIN) < 0.0) {
874  "Error: variance less than 0.0");
875  }
876  if (varStats.getValue(afwMath::MIN) == 0.0) {
878  "Error: variance equals 0.0, cannot inverse variance weight");
879  }
880 
881  lsst::afw::math::KernelList basisList =
882  boost::dynamic_pointer_cast<afwMath::LinearCombinationKernel>(this->_kernel)->getKernelList();
883 
884  unsigned int const nKernelParameters = basisList.size();
885  unsigned int const nBackgroundParameters = this->_fitForBackground ? 1 : 0;
886  unsigned int const nParameters = nKernelParameters + nBackgroundParameters;
887 
888  std::vector<boost::shared_ptr<afwMath::Kernel> >::const_iterator kiter = basisList.begin();
889 
890  /*
891  NOTE : If we are using these views in Afw's Image space, we need to
892  make sure and compensate for the XY0 of the image:
893 
894  afwGeom::Box2I fullBBox = templateImage.getBBox();
895  int maskStartCol = maskBox.getMinX();
896  int maskStartRow = maskBox.getMinY();
897  int maskEndCol = maskBox.getMaxX();
898  int maskEndRow = maskBox.getMaxY();
899 
900 
901  If we are going to be doing the slicing in Eigen matrices derived from
902  the images, ignore the XY0.
903 
904  afwGeom::Box2I fullBBox = templateImage.getBBox(afwImage::LOCAL);
905 
906  int maskStartCol = maskBox.getMinX() - templateImage.getX0();
907  int maskStartRow = maskBox.getMinY() - templateImage.getY0();
908  int maskEndCol = maskBox.getMaxX() - templateImage.getX0();
909  int maskEndRow = maskBox.getMaxY() - templateImage.getY0();
910 
911  */
912 
913 
914  afwGeom::Box2I shrunkBBox = (*kiter)->shrinkBBox(templateImage.getBBox());
915 
916  pexLog::TTrace<5>("lsst.ip.diffim.MaskedKernelSolution.build",
917  "Limits of good pixels after convolution: %d,%d -> %d,%d",
918  shrunkBBox.getMinX(), shrunkBBox.getMinY(),
919  shrunkBBox.getMaxX(), shrunkBBox.getMaxY());
920 
921  unsigned int const startCol = shrunkBBox.getMinX();
922  unsigned int const startRow = shrunkBBox.getMinY();
923  unsigned int const endCol = shrunkBBox.getMaxX();
924  unsigned int const endRow = shrunkBBox.getMaxY();
925 
926  /* NOTE: no endCol/endRow += 1 for block slicing, since we are doing the
927  slicing using Afw, not Eigen
928 
929  Eigen arrays have index 0,0 in the upper right, while LSST images
930  have 0,0 in the lower left. The y-axis is flipped. When translating
931  Images to Eigen matrices in ipDiffim::imageToEigenMatrix this is
932  accounted for. However, we still need to be aware of this fact if
933  addressing subregions of an Eigen matrix. This is why the slicing is
934  done in Afw, its cleaner.
935 
936  Please see examples/maskedKernel.cc for elaboration on some of the
937  tests done to make sure this indexing gets done right.
938 
939  */
940 
941 
942  /* Inner limits; of mask box */
943  int maskStartCol = maskBox.getMinX();
944  int maskStartRow = maskBox.getMinY();
945  int maskEndCol = maskBox.getMaxX();
946  int maskEndRow = maskBox.getMaxY();
947 
948  /*
949 
950  |---------------------------|
951  | Kernel Boundary |
952  | |---------------------| |
953  | | Top | |
954  | |......_________......| |
955  | | | | | |
956  | | L | Box | R | |
957  | | | | | |
958  | |......---------......| |
959  | | Bottom | |
960  | |---------------------| |
961  | |
962  |---------------------------|
963 
964  4 regions we want to extract from the pixels: top bottom left right
965 
966  */
967  afwGeom::Box2I tBox = afwGeom::Box2I(afwGeom::Point2I(startCol, maskEndRow + 1),
968  afwGeom::Point2I(endCol, endRow));
969 
970  afwGeom::Box2I bBox = afwGeom::Box2I(afwGeom::Point2I(startCol, startRow),
971  afwGeom::Point2I(endCol, maskStartRow - 1));
972 
973  afwGeom::Box2I lBox = afwGeom::Box2I(afwGeom::Point2I(startCol, maskStartRow),
974  afwGeom::Point2I(maskStartCol - 1, maskEndRow));
975 
976  afwGeom::Box2I rBox = afwGeom::Box2I(afwGeom::Point2I(maskEndCol + 1, maskStartRow),
977  afwGeom::Point2I(endCol, maskEndRow));
978 
979  pexLog::TTrace<5>("lsst.ip.diffim.MaskedKernelSolution.build",
980  "Upper good pixel region: %d,%d -> %d,%d",
981  tBox.getMinX(), tBox.getMinY(), tBox.getMaxX(), tBox.getMaxY());
982  pexLog::TTrace<5>("lsst.ip.diffim.MaskedKernelSolution.build",
983  "Bottom good pixel region: %d,%d -> %d,%d",
984  bBox.getMinX(), bBox.getMinY(), bBox.getMaxX(), bBox.getMaxY());
985  pexLog::TTrace<5>("lsst.ip.diffim.MaskedKernelSolution.build",
986  "Left good pixel region: %d,%d -> %d,%d",
987  lBox.getMinX(), lBox.getMinY(), lBox.getMaxX(), lBox.getMaxY());
988  pexLog::TTrace<5>("lsst.ip.diffim.MaskedKernelSolution.build",
989  "Right good pixel region: %d,%d -> %d,%d",
990  rBox.getMinX(), rBox.getMinY(), rBox.getMaxX(), rBox.getMaxY());
991 
992  std::vector<afwGeom::Box2I> boxArray;
993  boxArray.push_back(tBox);
994  boxArray.push_back(bBox);
995  boxArray.push_back(lBox);
996  boxArray.push_back(rBox);
997 
998  int totalSize = tBox.getWidth() * tBox.getHeight();
999  totalSize += bBox.getWidth() * bBox.getHeight();
1000  totalSize += lBox.getWidth() * lBox.getHeight();
1001  totalSize += rBox.getWidth() * rBox.getHeight();
1002 
1003  Eigen::MatrixXd eigenTemplate(totalSize, 1);
1004  Eigen::MatrixXd eigenScience(totalSize, 1);
1005  Eigen::MatrixXd eigeniVariance(totalSize, 1);
1006  eigenTemplate.setZero();
1007  eigenScience.setZero();
1008  eigeniVariance.setZero();
1009 
1010  boost::timer t;
1011  t.restart();
1012 
1013  int nTerms = 0;
1014  typename std::vector<afwGeom::Box2I>::iterator biter = boxArray.begin();
1015  for (; biter != boxArray.end(); ++biter) {
1016  int area = (*biter).getWidth() * (*biter).getHeight();
1017 
1018  afwImage::Image<InputT> siTemplate(templateImage, *biter);
1019  afwImage::Image<InputT> siScience(scienceImage, *biter);
1020  afwImage::Image<InputT> sVarEstimate(varianceEstimate, *biter);
1021 
1022  Eigen::MatrixXd eTemplate = imageToEigenMatrix(siTemplate);
1023  Eigen::MatrixXd eScience = imageToEigenMatrix(siScience);
1024  Eigen::MatrixXd eiVarEstimate = imageToEigenMatrix(sVarEstimate).array().inverse().matrix();
1025 
1026  eTemplate.resize(area, 1);
1027  eScience.resize(area, 1);
1028  eiVarEstimate.resize(area, 1);
1029 
1030  eigenTemplate.block(nTerms, 0, area, 1) = eTemplate.block(0, 0, area, 1);
1031  eigenScience.block(nTerms, 0, area, 1) = eScience.block(0, 0, area, 1);
1032  eigeniVariance.block(nTerms, 0, area, 1) = eiVarEstimate.block(0, 0, area, 1);
1033 
1034  nTerms += area;
1035  }
1036 
1037  afwImage::Image<InputT> cimage(templateImage.getDimensions());
1038 
1039  std::vector<boost::shared_ptr<Eigen::MatrixXd> > convolvedEigenList(nKernelParameters);
1040  typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiter =
1041  convolvedEigenList.begin();
1042  /* Create C_i in the formalism of Alard & Lupton */
1043  for (kiter = basisList.begin(); kiter != basisList.end(); ++kiter, ++eiter) {
1044  afwMath::convolve(cimage, templateImage, **kiter, false); /* cimage stores convolved image */
1045  Eigen::MatrixXd cMat(totalSize, 1);
1046  cMat.setZero();
1047 
1048  int nTerms = 0;
1049  typename std::vector<afwGeom::Box2I>::iterator biter = boxArray.begin();
1050  for (; biter != boxArray.end(); ++biter) {
1051  int area = (*biter).getWidth() * (*biter).getHeight();
1052 
1053  afwImage::Image<InputT> csubimage(cimage, *biter);
1054  Eigen::MatrixXd esubimage = imageToEigenMatrix(csubimage);
1055  esubimage.resize(area, 1);
1056  cMat.block(nTerms, 0, area, 1) = esubimage.block(0, 0, area, 1);
1057 
1058  nTerms += area;
1059  }
1060 
1061  boost::shared_ptr<Eigen::MatrixXd> cMatPtr (new Eigen::MatrixXd(cMat));
1062  *eiter = cMatPtr;
1063 
1064  }
1065 
1066  double time = t.elapsed();
1067  pexLog::TTrace<5>("lsst.ip.diffim.MaskedKernelSolution.build",
1068  "Total compute time to do basis convolutions : %.2f s", time);
1069  t.restart();
1070 
1071  /*
1072  Load matrix with all values from convolvedEigenList : all images
1073  (eigeniVariance, convolvedEigenList) must be the same size
1074  */
1075  Eigen::MatrixXd cMat(eigenTemplate.col(0).size(), nParameters);
1076  typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiterj =
1077  convolvedEigenList.begin();
1078  typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiterE =
1079  convolvedEigenList.end();
1080  for (unsigned int kidxj = 0; eiterj != eiterE; eiterj++, kidxj++) {
1081  cMat.col(kidxj) = (*eiterj)->col(0);
1082  }
1083  /* Treat the last "image" as all 1's to do the background calculation. */
1084  if (this->_fitForBackground)
1085  cMat.col(nParameters-1).fill(1.);
1086 
1087  this->_cMat.reset(new Eigen::MatrixXd(cMat));
1088  this->_ivVec.reset(new Eigen::VectorXd(eigeniVariance.col(0)));
1089  this->_iVec.reset(new Eigen::VectorXd(eigenScience.col(0)));
1090 
1091  /* Make these outside of solve() so I can check condition number */
1092  this->_mMat.reset(
1093  new Eigen::MatrixXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_cMat))
1094  );
1095  this->_bVec.reset(
1096  new Eigen::VectorXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_iVec))
1097  );
1098 
1099  }
1100  /*******************************************************************************************************/
1101 
1102 
1103  template <typename InputT>
1105  lsst::afw::math::KernelList const& basisList,
1106  bool fitForBackground,
1107  boost::shared_ptr<Eigen::MatrixXd> hMat,
1109  )
1110  :
1111  StaticKernelSolution<InputT>(basisList, fitForBackground),
1112  _hMat(hMat),
1113  _policy(policy)
1114  {};
1115 
1116  template <typename InputT>
1118  Eigen::MatrixXd vMat = (this->_cMat)->jacobiSvd().matrixV();
1119  Eigen::MatrixXd vMatvMatT = vMat * vMat.transpose();
1120 
1121  /* Find pseudo inverse of mMat, which may be ill conditioned */
1122  Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(*(this->_mMat));
1123  Eigen::MatrixXd const& rMat = eVecValues.eigenvectors();
1124  Eigen::VectorXd eValues = eVecValues.eigenvalues();
1125  double eMax = eValues.maxCoeff();
1126  for (int i = 0; i != eValues.rows(); ++i) {
1127  if (eValues(i) == 0.0) {
1128  eValues(i) = 0.0;
1129  }
1130  else if ((eMax / eValues(i)) > maxCond) {
1131  pexLog::TTrace<5>("lsst.ip.diffim.RegularizedKernelSolution.estimateRisk",
1132  "Truncating eValue %d; %.5e / %.5e = %.5e vs. %.5e",
1133  i, eMax, eValues(i), eMax / eValues(i), maxCond);
1134  eValues(i) = 0.0;
1135  }
1136  else {
1137  eValues(i) = 1.0 / eValues(i);
1138  }
1139  }
1140  Eigen::MatrixXd mInv = rMat * eValues.asDiagonal() * rMat.transpose();
1141 
1142  std::vector<double> lambdas = _createLambdaSteps();
1143  std::vector<double> risks;
1144  for (unsigned int i = 0; i < lambdas.size(); i++) {
1145  double l = lambdas[i];
1146  Eigen::MatrixXd mLambda = *(this->_mMat) + l * (*_hMat);
1147 
1148  try {
1149  KernelSolution::solve(mLambda, *(this->_bVec));
1150  } catch (pexExcept::Exception &e) {
1151  LSST_EXCEPT_ADD(e, "Unable to solve regularized kernel matrix");
1152  throw e;
1153  }
1154  Eigen::VectorXd term1 = (this->_aVec->transpose() * vMatvMatT * *(this->_aVec));
1155  if (term1.size() != 1)
1156  throw LSST_EXCEPT(pexExcept::Exception, "Matrix size mismatch");
1157 
1158  double term2a = (vMatvMatT * mLambda.inverse()).trace();
1159 
1160  Eigen::VectorXd term2b = (this->_aVec->transpose() * (mInv * *(this->_bVec)));
1161  if (term2b.size() != 1)
1162  throw LSST_EXCEPT(pexExcept::Exception, "Matrix size mismatch");
1163 
1164  double risk = term1(0) + 2 * (term2a - term2b(0));
1165  pexLog::TTrace<6>("lsst.ip.diffim.RegularizedKernelSolution.estimateRisk",
1166  "Lambda = %.3f, Risk = %.5e",
1167  l, risk);
1168  pexLog::TTrace<7>("lsst.ip.diffim.RegularizedKernelSolution.estimateRisk",
1169  "%.5e + 2 * (%.5e - %.5e)",
1170  term1(0), term2a, term2b(0));
1171  risks.push_back(risk);
1172  }
1173  std::vector<double>::iterator it = min_element(risks.begin(), risks.end());
1174  int index = distance(risks.begin(), it);
1175  pexLog::TTrace<5>("lsst.ip.diffim.RegularizedKernelSolution.estimateRisk",
1176  "Minimum Risk = %.3e at lambda = %.3e", risks[index], lambdas[index]);
1177 
1178  return lambdas[index];
1179 
1180  }
1181 
1182  template <typename InputT>
1183  boost::shared_ptr<Eigen::MatrixXd> RegularizedKernelSolution<InputT>::getM(bool includeHmat) {
1184  if (includeHmat == true) {
1185  return (boost::shared_ptr<Eigen::MatrixXd>(
1186  new Eigen::MatrixXd(*(this->_mMat) + _lambda * (*_hMat))
1187  ));
1188  }
1189  else {
1190  return this->_mMat;
1191  }
1192  }
1193 
1194  template <typename InputT>
1196 
1197  pexLog::TTrace<5>("lsst.ip.diffim.RegularizedKernelSolution.solve",
1198  "cMat is %d x %d; vVec is %d; iVec is %d; hMat is %d x %d",
1199  (*this->_cMat).rows(), (*this->_cMat).cols(), (*this->_ivVec).size(),
1200  (*this->_iVec).size(), (*_hMat).rows(), (*_hMat).cols());
1201 
1202  if (DEBUG_MATRIX2) {
1203  std::cout << "ID: " << (this->_id) << std::endl;
1204  std::cout << "C:" << std::endl;
1205  std::cout << (*this->_cMat) << std::endl;
1206  std::cout << "Sigma^{-1}:" << std::endl;
1207  std::cout << Eigen::MatrixXd((*this->_ivVec).asDiagonal()) << std::endl;
1208  std::cout << "Y:" << std::endl;
1209  std::cout << (*this->_iVec) << std::endl;
1210  std::cout << "H:" << std::endl;
1211  std::cout << (*_hMat) << std::endl;
1212  }
1213 
1214 
1215  this->_mMat.reset(
1216  new Eigen::MatrixXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_cMat))
1217  );
1218  this->_bVec.reset(
1219  new Eigen::VectorXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_iVec))
1220  );
1221 
1222 
1223  /* See N.R. 18.5
1224 
1225  Matrix equation to solve is Y = C a + N
1226  Y = vectorized version of I (I = image to not convolve)
1227  C_i = K_i (x) R (R = image to convolve)
1228  a = kernel coefficients
1229  N = noise
1230 
1231  If we reweight everything by the inverse square root of the noise
1232  covariance, we get a linear model with the identity matrix for
1233  the noise. The problem can then be solved using least squares,
1234  with the normal equations
1235 
1236  C^T Y = C^T C a
1237 
1238  or
1239 
1240  b = M a
1241 
1242  with
1243 
1244  b = C^T Y
1245  M = C^T C
1246  a = (C^T C)^{-1} C^T Y
1247 
1248 
1249  We can regularize the least square problem
1250 
1251  Y = C a + lambda a^T H a (+ N, which can be ignored)
1252 
1253  or the normal equations
1254 
1255  (C^T C + lambda H) a = C^T Y
1256 
1257 
1258  The solution to the regularization of the least squares problem is
1259 
1260  a = (C^T C + lambda H)^{-1} C^T Y
1261 
1262  The approximation to Y is
1263 
1264  C a = C (C^T C + lambda H)^{-1} C^T Y
1265 
1266  with smoothing matrix
1267 
1268  S = C (C^T C + lambda H)^{-1} C^T
1269 
1270  */
1271 
1272  std::string lambdaType = _policy.getString("lambdaType");
1273  if (lambdaType == "absolute") {
1274  _lambda = _policy.getDouble("lambdaValue");
1275  }
1276  else if (lambdaType == "relative") {
1277  _lambda = this->_mMat->trace() / this->_hMat->trace();
1278  _lambda *= _policy.getDouble("lambdaScaling");
1279  }
1280  else if (lambdaType == "minimizeBiasedRisk") {
1281  double tol = _policy.getDouble("maxConditionNumber");
1282  _lambda = estimateRisk(tol);
1283  }
1284  else if (lambdaType == "minimizeUnbiasedRisk") {
1285  _lambda = estimateRisk(std::numeric_limits<double>::max());
1286  }
1287  else {
1288  throw LSST_EXCEPT(pexExcept::Exception, "lambdaType in Policy not recognized");
1289  }
1290 
1291  pexLog::TTrace<5>("lsst.ip.diffim.KernelCandidate.build",
1292  "Applying kernel regularization with lambda = %.2e", _lambda);
1293 
1294 
1295  try {
1296  KernelSolution::solve(*(this->_mMat) + _lambda * (*_hMat), *(this->_bVec));
1297  } catch (pexExcept::Exception &e) {
1298  LSST_EXCEPT_ADD(e, "Unable to solve static kernel matrix");
1299  throw e;
1300  }
1301  /* Turn matrices into _kernel and _background */
1303  }
1304 
1305  template <typename InputT>
1307  std::vector<double> lambdas;
1308 
1309  std::string lambdaStepType = _policy.getString("lambdaStepType");
1310  if (lambdaStepType == "linear") {
1311  double lambdaLinMin = _policy.getDouble("lambdaLinMin");
1312  double lambdaLinMax = _policy.getDouble("lambdaLinMax");
1313  double lambdaLinStep = _policy.getDouble("lambdaLinStep");
1314  for (double l = lambdaLinMin; l <= lambdaLinMax; l += lambdaLinStep) {
1315  lambdas.push_back(l);
1316  }
1317  }
1318  else if (lambdaStepType == "log") {
1319  double lambdaLogMin = _policy.getDouble("lambdaLogMin");
1320  double lambdaLogMax = _policy.getDouble("lambdaLogMax");
1321  double lambdaLogStep = _policy.getDouble("lambdaLogStep");
1322  for (double l = lambdaLogMin; l <= lambdaLogMax; l += lambdaLogStep) {
1323  lambdas.push_back(pow(10, l));
1324  }
1325  }
1326  else {
1327  throw LSST_EXCEPT(pexExcept::Exception, "lambdaStepType in Policy not recognized");
1328  }
1329  return lambdas;
1330  }
1331 
1332  /*******************************************************************************************************/
1333 
1335  lsst::afw::math::KernelList const& basisList,
1336  lsst::afw::math::Kernel::SpatialFunctionPtr spatialKernelFunction,
1339  ) :
1340  KernelSolution(),
1341  _spatialKernelFunction(spatialKernelFunction),
1342  _constantFirstTerm(false),
1343  _kernel(),
1344  _background(background),
1345  _kSum(0.0),
1346  _policy(policy),
1347  _nbases(0),
1348  _nkt(0),
1349  _nbt(0),
1350  _nt(0) {
1351 
1352  bool isAlardLupton = _policy.getString("kernelBasisSet") == "alard-lupton";
1353  bool usePca = _policy.getBool("usePcaForSpatialKernel");
1354  if (isAlardLupton || usePca) {
1355  _constantFirstTerm = true;
1356  }
1357  this->_fitForBackground = _policy.getBool("fitForBackground");
1358 
1359  _nbases = basisList.size();
1360  _nkt = _spatialKernelFunction->getParameters().size();
1361  _nbt = _fitForBackground ? _background->getParameters().size() : 0;
1362  _nt = 0;
1363  if (_constantFirstTerm) {
1364  _nt = (_nbases - 1) * _nkt + 1 + _nbt;
1365  } else {
1366  _nt = _nbases * _nkt + _nbt;
1367  }
1368 
1369  boost::shared_ptr<Eigen::MatrixXd> mMat (new Eigen::MatrixXd(_nt, _nt));
1370  boost::shared_ptr<Eigen::VectorXd> bVec (new Eigen::VectorXd(_nt));
1371  (*mMat).setZero();
1372  (*bVec).setZero();
1373 
1374  this->_mMat = mMat;
1375  this->_bVec = bVec;
1376 
1379  );
1380 
1381  pexLog::TTrace<5>("lsst.ip.diffim.SpatialKernelSolution",
1382  "Initializing with size %d %d %d and constant first term = %s",
1383  _nkt, _nbt, _nt,
1384  _constantFirstTerm ? "true" : "false");
1385 
1386  }
1387 
1388  void SpatialKernelSolution::addConstraint(float xCenter, float yCenter,
1389  boost::shared_ptr<Eigen::MatrixXd> qMat,
1390  boost::shared_ptr<Eigen::VectorXd> wVec) {
1391 
1392  pexLog::TTrace<8>("lsst.ip.diffim.SpatialKernelSolution.addConstraint",
1393  "Adding candidate at %f, %f", xCenter, yCenter);
1394 
1395  /* Calculate P matrices */
1396  /* Pure kernel terms */
1397  Eigen::VectorXd pK(_nkt);
1398  std::vector<double> paramsK = _spatialKernelFunction->getParameters();
1399  for (int idx = 0; idx < _nkt; idx++) { paramsK[idx] = 0.0; }
1400  for (int idx = 0; idx < _nkt; idx++) {
1401  paramsK[idx] = 1.0;
1402  _spatialKernelFunction->setParameters(paramsK);
1403  pK(idx) = (*_spatialKernelFunction)(xCenter, yCenter); /* Assume things don't vary over stamp */
1404  paramsK[idx] = 0.0;
1405  }
1406  Eigen::MatrixXd pKpKt = (pK * pK.transpose());
1407 
1408  Eigen::VectorXd pB;
1409  Eigen::MatrixXd pBpBt;
1410  Eigen::MatrixXd pKpBt;
1411  if (_fitForBackground) {
1412  pB = Eigen::VectorXd(_nbt);
1413 
1414  /* Pure background terms */
1415  std::vector<double> paramsB = _background->getParameters();
1416  for (int idx = 0; idx < _nbt; idx++) { paramsB[idx] = 0.0; }
1417  for (int idx = 0; idx < _nbt; idx++) {
1418  paramsB[idx] = 1.0;
1419  _background->setParameters(paramsB);
1420  pB(idx) = (*_background)(xCenter, yCenter); /* Assume things don't vary over stamp */
1421  paramsB[idx] = 0.0;
1422  }
1423  pBpBt = (pB * pB.transpose());
1424 
1425  /* Cross terms */
1426  pKpBt = (pK * pB.transpose());
1427  }
1428 
1429  if (DEBUG_MATRIX) {
1430  std::cout << "Spatial weights" << std::endl;
1431  std::cout << "pKpKt " << pKpKt << std::endl;
1432  if (_fitForBackground) {
1433  std::cout << "pBpBt " << pBpBt << std::endl;
1434  std::cout << "pKpBt " << pKpBt << std::endl;
1435  }
1436  }
1437 
1438  if (DEBUG_MATRIX) {
1439  std::cout << "Spatial matrix inputs" << std::endl;
1440  std::cout << "M " << (*qMat) << std::endl;
1441  std::cout << "B " << (*wVec) << std::endl;
1442  }
1443 
1444  /* first index to start the spatial blocks; default=0 for non-constant first term */
1445  int m0 = 0;
1446  /* how many rows/cols to adjust the matrices/vectors; default=0 for non-constant first term */
1447  int dm = 0;
1448  /* where to start the background terms; this is always true */
1449  int mb = _nt - _nbt;
1450 
1451  if (_constantFirstTerm) {
1452  m0 = 1; /* we need to manually fill in the first (non-spatial) terms below */
1453  dm = _nkt-1; /* need to shift terms due to lack of spatial variation in first term */
1454 
1455  (*_mMat)(0, 0) += (*qMat)(0,0);
1456  for(int m2 = 1; m2 < _nbases; m2++) {
1457  (*_mMat).block(0, m2*_nkt-dm, 1, _nkt) += (*qMat)(0,m2) * pK.transpose();
1458  }
1459  (*_bVec)(0) += (*wVec)(0);
1460 
1461  if (_fitForBackground) {
1462  (*_mMat).block(0, mb, 1, _nbt) += (*qMat)(0,_nbases) * pB.transpose();
1463  }
1464  }
1465 
1466  /* Fill in the spatial blocks */
1467  for(int m1 = m0; m1 < _nbases; m1++) {
1468  /* Diagonal kernel-kernel term; only use upper triangular part of pKpKt */
1469  (*_mMat).block(m1*_nkt-dm, m1*_nkt-dm, _nkt, _nkt) +=
1470  (pKpKt * (*qMat)(m1,m1)).triangularView<Eigen::Upper>();
1471 
1472  /* Kernel-kernel terms */
1473  for(int m2 = m1+1; m2 < _nbases; m2++) {
1474  (*_mMat).block(m1*_nkt-dm, m2*_nkt-dm, _nkt, _nkt) += (*qMat)(m1,m2) * pKpKt;
1475  }
1476 
1477  if (_fitForBackground) {
1478  /* Kernel cross terms with background */
1479  (*_mMat).block(m1*_nkt-dm, mb, _nkt, _nbt) += (*qMat)(m1,_nbases) * pKpBt;
1480  }
1481 
1482  /* B vector */
1483  (*_bVec).segment(m1*_nkt-dm, _nkt) += (*wVec)(m1) * pK;
1484  }
1485 
1486  if (_fitForBackground) {
1487  /* Background-background terms only */
1488  (*_mMat).block(mb, mb, _nbt, _nbt) +=
1489  (pBpBt * (*qMat)(_nbases,_nbases)).triangularView<Eigen::Upper>();
1490  (*_bVec).segment(mb, _nbt) += (*wVec)(_nbases) * pB;
1491  }
1492 
1493  if (DEBUG_MATRIX) {
1494  std::cout << "Spatial matrix outputs" << std::endl;
1495  std::cout << "mMat " << (*_mMat) << std::endl;
1496  std::cout << "bVec " << (*_bVec) << std::endl;
1497  }
1498 
1499  }
1500 
1503  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return image");
1504  }
1506  new afwImage::Image<afwMath::Kernel::Pixel>(_kernel->getDimensions())
1507  );
1508  (void)_kernel->computeImage(*image, false, pos[0], pos[1]);
1509  return image;
1510  }
1511 
1513  /* Fill in the other half of mMat */
1514  for (int i = 0; i < _nt; i++) {
1515  for (int j = i+1; j < _nt; j++) {
1516  (*_mMat)(j,i) = (*_mMat)(i,j);
1517  }
1518  }
1519 
1520  try {
1522  } catch (pexExcept::Exception &e) {
1523  LSST_EXCEPT_ADD(e, "Unable to solve spatial kernel matrix");
1524  throw e;
1525  }
1526  /* Turn matrices into _kernel and _background */
1527  _setKernel();
1528  }
1529 
1533  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return solution");
1534  }
1535 
1536  return std::make_pair(_kernel, _background);
1537  }
1538 
1540  /* Report on condition number */
1541  double cNumber = this->getConditionNumber(EIGENVALUE);
1542 
1543  if (_nkt == 1) {
1544  /* Not spatially varying; this fork is a specialization for convolution speed--up */
1545 
1546  /* Set the basis coefficients */
1547  std::vector<double> kCoeffs(_nbases);
1548  for (int i = 0; i < _nbases; i++) {
1549  if (std::isnan((*_aVec)(i))) {
1550  throw LSST_EXCEPT(
1552  str(boost::format(
1553  "I. Unable to determine spatial kernel solution %d (nan). Condition number = %.3e") % i % cNumber));
1554  }
1555  kCoeffs[i] = (*_aVec)(i);
1556  }
1557  lsst::afw::math::KernelList basisList =
1558  boost::dynamic_pointer_cast<afwMath::LinearCombinationKernel>(_kernel)->getKernelList();
1559  _kernel.reset(
1560  new afwMath::LinearCombinationKernel(basisList, kCoeffs)
1561  );
1562  }
1563  else {
1564 
1565  /* Set the kernel coefficients */
1566  std::vector<std::vector<double> > kCoeffs;
1567  kCoeffs.reserve(_nbases);
1568  for (int i = 0, idx = 0; i < _nbases; i++) {
1569  kCoeffs.push_back(std::vector<double>(_nkt));
1570 
1571  /* Deal with the possibility the first term doesn't vary spatially */
1572  if ((i == 0) && (_constantFirstTerm)) {
1573  if (std::isnan((*_aVec)(idx))) {
1574  throw LSST_EXCEPT(
1576  str(boost::format(
1577  "II. Unable to determine spatial kernel solution %d (nan). Condition number = %.3e") % idx % cNumber));
1578  }
1579  kCoeffs[i][0] = (*_aVec)(idx++);
1580  }
1581  else {
1582  for (int j = 0; j < _nkt; j++) {
1583  if (std::isnan((*_aVec)(idx))) {
1584  throw LSST_EXCEPT(
1586  str(boost::format(
1587  "III. Unable to determine spatial kernel solution %d (nan). Condition number = %.3e") % idx % cNumber));
1588  }
1589  kCoeffs[i][j] = (*_aVec)(idx++);
1590  }
1591  }
1592  }
1593  _kernel->setSpatialParameters(kCoeffs);
1594  }
1595 
1596  /* Set kernel Sum */
1597  ImageT::Ptr image (new ImageT(_kernel->getDimensions()));
1598  _kSum = _kernel->computeImage(*image, false);
1599 
1600  /* Set the background coefficients */
1601  std::vector<double> bgCoeffs(_fitForBackground ? _nbt : 1);
1602  if (_fitForBackground) {
1603  for (int i = 0; i < _nbt; i++) {
1604  int idx = _nt - _nbt + i;
1605  if (std::isnan((*_aVec)(idx))) {
1607  str(boost::format(
1608  "Unable to determine spatial background solution %d (nan)") % (idx)));
1609  }
1610  bgCoeffs[i] = (*_aVec)(idx);
1611  }
1612  }
1613  else {
1614  bgCoeffs[0] = 0.;
1615  }
1616  _background->setParameters(bgCoeffs);
1617  }
1618 
1619 /***********************************************************************************************************/
1620 //
1621 // Explicit instantiations
1622 //
1623  typedef float InputT;
1624 
1625  template class StaticKernelSolution<InputT>;
1626  template class MaskedKernelSolution<InputT>;
1627  template class RegularizedKernelSolution<InputT>;
1628 
1629 }}} // end of namespace lsst::ip::diffim
An include file to include the public header files for lsst::afw::math.
void _setKernel()
Set kernel after solution.
iterator end() const
Return an STL compliant iterator to the end of the image.
Definition: Image.cc:331
void flattenArray(Footprint const &fp, ndarray::Array< T, N, C > const &src, ndarray::Array< U, N-1, D > const &dest, lsst::afw::geom::Point2I const &xy0=lsst::afw::geom::Point2I())
Flatten the first two dimensions of an array.
boost::shared_ptr< lsst::afw::math::Function2< double > > SpatialFunctionPtr
Definition: Kernel.h:143
int getMaxY() const
Definition: Box.h:129
An include file to include the header files for lsst::afw::geom.
geom::Extent2I const getDimensions() const
Return the Kernel&#39;s dimensions (width, height)
Definition: Kernel.h:226
Eigen matrix objects that present a view into an ndarray::Array.
Eigen::MatrixXi maskToEigenMatrix(lsst::afw::image::Mask< lsst::afw::image::MaskPixel > const &mask)
boost::shared_ptr< Footprint > Ptr
Definition: Footprint.h:67
void writeFits(std::string const &fileName, boost::shared_ptr< lsst::daf::base::PropertySet const > metadata=boost::shared_ptr< lsst::daf::base::PropertySet >(), std::string const &mode="w") const
Write a mask to a regular FITS file.
estimate sample minimum
Definition: Statistics.h:76
Eigen3 view into an ndarray::Array.
Definition: eigen.h:232
virtual lsst::afw::image::Image< lsst::afw::math::Kernel::Pixel >::Ptr makeKernelImage()
boost::shared_ptr< FootprintList > getFootprints()
Definition: FootprintSet.h:146
bool getBool(const std::string &name) const
Definition: Policy.h:589
boost::uint16_t MaskPixel
a container for holding hierarchical configuration data in memory.
Definition: Policy.h:169
void setKernelParameters(std::vector< double > const &params)
Set the kernel parameters of a spatially invariant kernel.
Definition: Kernel.h:412
const std::string getString(const std::string &name) const
Definition: Policy.h:631
definition of the Trace messaging facilities
SpatialKernelSolution(lsst::afw::math::KernelList const &basisList, lsst::afw::math::Kernel::SpatialFunctionPtr spatialKernelFunction, lsst::afw::math::Kernel::SpatialFunctionPtr background, lsst::pex::policy::Policy policy)
boost::shared_ptr< LinearCombinationKernel > Ptr
Definition: Kernel.h:818
lsst::afw::math::Kernel::Ptr _kernel
Derived single-object convolution kernel.
Use (pixels &amp; (given mask))
Definition: Threshold.h:49
A Threshold is used to pass a threshold value to detection algorithms.
Definition: Threshold.h:44
virtual lsst::afw::math::Kernel::Ptr getKernel()
geom::Box2I getBBox(ImageOrigin origin=PARENT) const
Definition: Image.h:377
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
StaticKernelSolution(lsst::afw::math::KernelList const &basisList, bool fitForBackground)
int isnan(T t)
Definition: ieee.h:110
def trace
Definition: log.py:96
boost::shared_ptr< Eigen::MatrixXd > getM()
An integer coordinate rectangle.
Definition: Box.h:53
Declaration of classes to store the solution for convolution kernels.
table::Key< table::Array< Kernel::Pixel > > image
Definition: FixedKernel.cc:117
void addConstraint(float xCenter, float yCenter, boost::shared_ptr< Eigen::MatrixXd > qMat, boost::shared_ptr< Eigen::VectorXd > wVec)
double getValue(Property const prop=NOTHING) const
Return the value of the desired property (if specified in the constructor)
Definition: Statistics.cc:1009
An include file to include the header files for lsst::afw::image.
geom::Point2I getXY0() const
Definition: Image.h:264
boost::shared_ptr< Kernel > Ptr
Definition: Kernel.h:141
virtual void build(lsst::afw::image::Image< InputT > const &templateImage, lsst::afw::image::Image< InputT > const &scienceImage, lsst::afw::image::Image< lsst::afw::image::VariancePixel > const &varianceEstimate)
void _setKernel()
Set kernel after solution.
double _kSum
Derived kernel sum.
static int _SolutionId
Unique identifier for solution.
void _setKernelUncertainty()
Not implemented.
bool _fitForBackground
Background terms included in fit.
geom::Extent2I getDimensions() const
Return the image&#39;s size; useful for passing to constructors.
Definition: Image.h:298
An include file to include the header files for lsst::afw::detection.
lsst::afw::image::Image< lsst::afw::math::Kernel::Pixel > ImageT
int getMinY() const
Definition: Box.h:125
afwMath::LinearCombinationKernel const & _kernel
virtual double getConditionNumber(ConditionNumberType conditionType)
int _id
Definition: Mask.cc:679
lsst::afw::image::Image< lsst::afw::math::Kernel::Pixel >::Ptr makeKernelImage(lsst::afw::geom::Point2D const &pos)
double computeImage(lsst::afw::image::Image< Pixel > &image, bool doNormalize, double x=0.0, double y=0.0) const
Compute an image (pixellized representation of the kernel) in place.
Definition: Kernel.cc:94
A kernel that is a linear combination of fixed basis kernels.
Definition: Kernel.h:814
Vector< T, N > makeVector(T v1, T v2,..., T vN)
Variadic constructor for Vector.
int getMinX() const
Definition: Box.h:124
boost::shared_ptr< Image< PixelT > > Ptr
Definition: Image.h:418
boost::shared_ptr< Eigen::VectorXd > _bVec
Derived least squares B vector.
virtual void buildOrig(lsst::afw::image::Image< InputT > const &templateImage, lsst::afw::image::Image< InputT > const &scienceImage, lsst::afw::image::Image< lsst::afw::image::VariancePixel > const &varianceEstimate, lsst::afw::image::Mask< lsst::afw::image::MaskPixel > pixelMask)
lsst::afw::math::Kernel::SpatialFunctionPtr _spatialKernelFunction
Spatial function for Kernel.
static MaskPixelT getPlaneBitMask(const std::vector< std::string > &names)
Return the bitmask corresponding to a vector of plane names OR&#39;d together.
Definition: Mask.cc:860
A set of pixels in an Image.
Definition: Footprint.h:62
int _nbt
Number of background terms.
bool _constantFirstTerm
Is the first term constant.
lsst::pex::policy::Policy _policy
Policy to control processing.
int getWidth() const
Definition: Box.h:154
detail::SimpleInitializer< N > allocate(Vector< int, N > const &shape)
Create an expression that allocates uninitialized memory for an array.
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46
A set of Footprints, associated with a MaskedImage.
Definition: FootprintSet.h:53
lsst::afw::math::LinearCombinationKernel::Ptr _kernel
Spatial convolution kernel.
boost::shared_ptr< Eigen::VectorXd > _aVec
Derived least squares solution matrix.
MaskedKernelSolution(lsst::afw::math::KernelList const &basisList, bool fitForBackground)
int getHeight() const
Definition: Box.h:155
lsst::afw::math::Kernel::SpatialFunctionPtr _background
Spatial background model.
void convolve(OutImageT &convolvedImage, InImageT const &inImage, KernelT const &kernel, ConvolutionControl const &convolutionControl=ConvolutionControl())
Convolve an Image or MaskedImage with a Kernel, setting pixels of an existing output image...
RegularizedKernelSolution(lsst::afw::math::KernelList const &basisList, bool fitForBackground, boost::shared_ptr< Eigen::MatrixXd > hMat, lsst::pex::policy::Policy policy)
int _nkt
Number of kernel terms.
boost::shared_ptr< Eigen::MatrixXd > _mMat
Derived least squares M matrix.
std::pair< lsst::afw::math::LinearCombinationKernel::Ptr, lsst::afw::math::Kernel::SpatialFunctionPtr > getSolutionPair()
MaskT setMaskFromFootprint(lsst::afw::image::Mask< MaskT > *mask, Footprint const &footprint, MaskT const bitmask)
OR bitmask into all the Mask&#39;s pixels that are in the Footprint.
int _nbases
Number of basis functions.
Statistics makeStatistics(afwImage::Mask< afwImage::MaskPixel > const &msk, int const flags, StatisticsControl const &sctrl)
Specialization to handle Masks.
Definition: Statistics.cc:1082
virtual void buildSingleMaskOrig(lsst::afw::image::Image< InputT > const &templateImage, lsst::afw::image::Image< InputT > const &scienceImage, lsst::afw::image::Image< lsst::afw::image::VariancePixel > const &varianceEstimate, lsst::afw::geom::Box2I maskBox)
virtual std::pair< boost::shared_ptr< lsst::afw::math::Kernel >, double > getSolutionPair()
virtual void buildWithMask(lsst::afw::image::Image< InputT > const &templateImage, lsst::afw::image::Image< InputT > const &scienceImage, lsst::afw::image::Image< lsst::afw::image::VariancePixel > const &varianceEstimate, lsst::afw::image::Mask< lsst::afw::image::MaskPixel > const &pixelMask)
int getMaxX() const
Definition: Box.h:128
iterator begin() const
Definition: Image.cc:325
#define LSST_EXCEPT_ADD(e, m)
Definition: Exception.h:51
A class to represent a 2-dimensional array of pixels.
Definition: Image.h:415
#define DEBUG_MATRIX
#define DEBUG_MATRIX2
Image Subtraction helper functions.
MaskT setMaskFromFootprintList(lsst::afw::image::Mask< MaskT > *mask, std::vector< boost::shared_ptr< Footprint >> const &footprints, MaskT const bitmask)
OR bitmask into all the Mask&#39;s pixels which are in the set of Footprints.
Eigen::MatrixXd imageToEigenMatrix(lsst::afw::image::Image< PixelT > const &img)
Turns Image into a 2-D Eigen Matrix.
std::vector< boost::shared_ptr< Kernel > > KernelList
Definition: Kernel.h:542