16 #include "boost/shared_ptr.hpp"
17 #include "boost/timer.hpp"
20 #include "Eigen/Cholesky"
23 #include "Eigen/Eigenvalues"
30 #include "lsst/afw/detection/FootprintArray.cc"
40 #define DEBUG_MATRIX 0
41 #define DEBUG_MATRIX2 0
43 namespace afwDet = lsst::afw::detection;
44 namespace afwMath = lsst::afw::math;
45 namespace afwGeom = lsst::afw::geom;
47 namespace pexLog = lsst::pex::logging;
48 namespace pexExcept = lsst::pex::exceptions;
58 boost::shared_ptr<Eigen::MatrixXd> mMat,
59 boost::shared_ptr<Eigen::VectorXd> bVec,
67 _fitForBackground(fitForBackground)
78 _fitForBackground(fitForBackground)
87 _fitForBackground(true)
100 switch (conditionType) {
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);
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);
124 throw LSST_EXCEPT(pexExcept::InvalidParameterError,
125 "Undefined ConditionNumberType : only EIGENVALUE, SVD allowed.");
132 Eigen::VectorXd bVec) {
135 std::cout <<
"M " << std::endl;
136 std::cout << mMat << std::endl;
137 std::cout <<
"B " << std::endl;
138 std::cout << bVec << std::endl;
141 Eigen::VectorXd aVec = Eigen::VectorXd::Zero(bVec.size());
146 pexLog::TTrace<4>(
"lsst.ip.difim.KernelSolution.solve",
147 "Solving for kernel");
149 Eigen::FullPivLU<Eigen::MatrixXd> lu(mMat);
150 if (lu.isInvertible()) {
151 aVec = lu.solve(bVec);
153 pexLog::TTrace<5>(
"lsst.ip.diffim.KernelSolution.solve",
154 "Unable to determine kernel via LU");
159 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(mMat);
160 Eigen::MatrixXd
const& rMat = eVecValues.eigenvectors();
161 Eigen::VectorXd eValues = eVecValues.eigenvalues();
163 for (
int i = 0; i != eValues.rows(); ++i) {
164 if (eValues(i) != 0.0) {
165 eValues(i) = 1.0/eValues(i);
169 aVec = rMat * eValues.asDiagonal() * rMat.transpose() * bVec;
173 pexLog::TTrace<5>(
"lsst.ip.diffim.KernelSolution.solve",
174 "Unable to determine kernel via eigen-values");
180 double time = t.elapsed();
181 pexLog::TTrace<5>(
"lsst.ip.diffim.KernelSolution.solve",
182 "Compute time for matrix math : %.2f s", time);
185 std::cout <<
"A " << std::endl;
186 std::cout << aVec << std::endl;
189 _aVec = boost::shared_ptr<Eigen::VectorXd>(
new Eigen::VectorXd(aVec));
194 template <
typename InputT>
197 bool fitForBackground
208 std::vector<double> kValues(basisList.size());
209 _kernel = boost::shared_ptr<afwMath::Kernel>(
214 template <
typename InputT>
222 template <
typename InputT>
234 template <
typename InputT>
242 template <
typename InputT>
250 template <
typename InputT>
251 std::pair<boost::shared_ptr<lsst::afw::math::Kernel>,
double>
257 return std::make_pair(
_kernel, _background);
260 template <
typename InputT>
270 "Error: variance less than 0.0");
274 "Error: variance equals 0.0, cannot inverse variance weight");
280 unsigned int const nKernelParameters = basisList.size();
281 unsigned int const nBackgroundParameters = _fitForBackground ? 1 : 0;
282 unsigned int const nParameters = nKernelParameters + nBackgroundParameters;
284 std::vector<boost::shared_ptr<afwMath::Kernel> >::const_iterator kiter = basisList.begin();
319 unsigned int const startCol = goodBBox.
getMinX();
320 unsigned int const startRow = goodBBox.
getMinY();
322 unsigned int endCol = goodBBox.
getMaxX() + 1;
323 unsigned int endRow = goodBBox.
getMaxY() + 1;
338 startRow, startCol, endRow-startRow, endCol-startCol
339 ).array().inverse().matrix();
342 eigenTemplate.resize(eigenTemplate.rows()*eigenTemplate.cols(), 1);
343 eigenScience.resize(eigenScience.rows()*eigenScience.cols(), 1);
344 eigeniVariance.resize(eigeniVariance.rows()*eigeniVariance.cols(), 1);
350 std::vector<boost::shared_ptr<Eigen::MatrixXd> > convolvedEigenList(nKernelParameters);
353 typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiter =
354 convolvedEigenList.
begin();
356 for (kiter = basisList.begin(); kiter != basisList.end(); ++kiter, ++eiter) {
359 boost::shared_ptr<Eigen::MatrixXd> cMat (
365 cMat->resize(cMat->rows()*cMat->cols(), 1);
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);
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);
388 if (_fitForBackground)
389 cMat.col(nParameters-1).fill(1.);
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)));
396 _mMat.reset(
new Eigen::MatrixXd((*_cMat).transpose() * ((*_ivVec).asDiagonal() * (*_cMat))));
397 _bVec.reset(
new Eigen::VectorXd((*_cMat).transpose() * ((*_ivVec).asDiagonal() * (*_iVec))));
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());
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;
432 template <
typename InputT>
438 unsigned int const nParameters = _aVec->size();
439 unsigned int const nBackgroundParameters = _fitForBackground ? 1 : 0;
440 unsigned int const nKernelParameters =
442 if (nParameters != (nKernelParameters + nBackgroundParameters))
446 std::vector<double> kValues(nKernelParameters);
447 for (
unsigned int idx = 0; idx < nKernelParameters; idx++) {
450 str(
boost::format(
"Unable to determine kernel solution %d (nan)") % idx));
452 kValues[idx] = (*_aVec)(idx);
461 if (_fitForBackground) {
464 str(
boost::format(
"Unable to determine background solution %d (nan)") %
467 _background = (*_aVec)(nParameters-1);
472 template <
typename InputT>
502 template <
typename InputT>
505 bool fitForBackground
510 template <
typename InputT>
521 "Error: variance less than 0.0");
525 "Error: variance equals 0.0, cannot inverse variance weight");
533 std::vector<boost::shared_ptr<afwMath::Kernel> >::const_iterator kiter = basisList.begin();
547 int growPix = (*kiter)->getCtr().getX();
551 for (
typename afwDet::FootprintSet::FootprintList::iterator
568 finalMask.writeFits(
"finalmask.fits");
574 maskArray, templateImage.
getXY0());
580 arrayTemplate, templateImage.
getXY0());
586 arrayScience, scienceImage.
getXY0());
592 arrayVariance, varianceEstimate.
getXY0());
596 for (
int i = 0; i < maskEigen.size(); i++) {
597 if (maskEigen(i) == 0.0)
601 Eigen::VectorXd eigenTemplate(nGood);
602 Eigen::VectorXd eigenScience(nGood);
603 Eigen::VectorXd eigenVariance(nGood);
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);
618 unsigned int const nKernelParameters = basisList.size();
619 unsigned int const nBackgroundParameters = this->_fitForBackground ? 1 : 0;
620 unsigned int const nParameters = nKernelParameters + nBackgroundParameters;
626 std::vector<boost::shared_ptr<Eigen::VectorXd> >
627 convolvedEigenList(nKernelParameters);
630 typename std::vector<boost::shared_ptr<Eigen::VectorXd> >::iterator eiter =
631 convolvedEigenList.
begin();
634 for (kiter = basisList.begin(); kiter != basisList.end(); ++kiter, ++eiter) {
640 arrayC, cimage.getXY0());
643 Eigen::VectorXd eigenC(nGood);
645 for (
int i = 0; i < maskEigen.size(); i++) {
646 if (maskEigen(i) == 0.0) {
647 eigenC(nUsed) = eigenC0(i);
652 boost::shared_ptr<Eigen::VectorXd>
653 eigenCptr (
new Eigen::VectorXd(eigenC));
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);
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);
673 if (this->_fitForBackground)
674 cMat.col(nParameters-1).fill(1.);
676 this->_cMat.reset(
new Eigen::MatrixXd(cMat));
679 this->_ivVec.reset(
new Eigen::VectorXd(eigenVariance.array().inverse().matrix()));
680 this->_iVec.reset(
new Eigen::VectorXd(eigenScience));
684 new Eigen::MatrixXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_cMat))
687 new Eigen::VectorXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_iVec))
692 template <
typename InputT>
703 "Error: variance less than 0.0");
707 "Error: variance equals 0.0, cannot inverse variance weight");
712 std::vector<boost::shared_ptr<afwMath::Kernel> >::const_iterator kiter = basisList.begin();
723 unsigned int const nKernelParameters = basisList.size();
724 unsigned int const nBackgroundParameters = this->_fitForBackground ? 1 : 0;
725 unsigned int const nParameters = nKernelParameters + nBackgroundParameters;
733 pexLog::TTrace<5>(
"lsst.ip.diffim.MaskedKernelSolution.build",
734 "Limits of good pixels after convolution: %d,%d -> %d,%d (local)",
739 unsigned int startCol = shrunkLocalBBox.
getMinX();
740 unsigned int startRow = shrunkLocalBBox.
getMinY();
741 unsigned int endCol = shrunkLocalBBox.
getMaxX();
742 unsigned int endRow = shrunkLocalBBox.
getMaxY();
765 startRow, startCol, endRow-startRow, endCol-startCol
766 ).array().inverse().matrix();
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);
775 Eigen::MatrixXd maskedEigenTemplate(eigenTemplate.rows(), 1);
776 Eigen::MatrixXd maskedEigenScience(eigenScience.rows(), 1);
777 Eigen::MatrixXd maskedEigeniVariance(eigeniVariance.rows(), 1);
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);
788 eigenTemplate = maskedEigenTemplate.block(0, 0, nGood, 1);
789 eigenScience = maskedEigenScience.block(0, 0, nGood, 1);
790 eigeniVariance = maskedEigeniVariance.block(0, 0, nGood, 1);
797 std::vector<boost::shared_ptr<Eigen::MatrixXd> > convolvedEigenList(nKernelParameters);
800 typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiter =
801 convolvedEigenList.
begin();
803 for (kiter = basisList.begin(); kiter != basisList.end(); ++kiter, ++eiter) {
810 cMat.resize(cMat.rows() * cMat.cols(), 1);
813 Eigen::MatrixXd maskedcMat(cMat.rows(), 1);
815 for (
int i = 0; i < eMask.rows(); i++) {
816 if (eMask(i, 0) == 0) {
817 maskedcMat(nGood, 0) = cMat(i, 0);
821 cMat = maskedcMat.block(0, 0, nGood, 1);
822 boost::shared_ptr<Eigen::MatrixXd> cMatPtr (
new Eigen::MatrixXd(cMat));
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);
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);
844 if (this->_fitForBackground)
845 cMat.col(nParameters-1).fill(1.);
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)));
853 new Eigen::MatrixXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_cMat))
856 new Eigen::VectorXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_iVec))
863 template <
typename InputT>
874 "Error: variance less than 0.0");
878 "Error: variance equals 0.0, cannot inverse variance weight");
884 unsigned int const nKernelParameters = basisList.size();
885 unsigned int const nBackgroundParameters = this->_fitForBackground ? 1 : 0;
886 unsigned int const nParameters = nKernelParameters + nBackgroundParameters;
888 std::vector<boost::shared_ptr<afwMath::Kernel> >::const_iterator kiter = basisList.begin();
916 pexLog::TTrace<5>(
"lsst.ip.diffim.MaskedKernelSolution.build",
917 "Limits of good pixels after convolution: %d,%d -> %d,%d",
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();
943 int maskStartCol = maskBox.
getMinX();
944 int maskStartRow = maskBox.
getMinY();
945 int maskEndCol = maskBox.
getMaxX();
946 int maskEndRow = maskBox.
getMaxY();
979 pexLog::TTrace<5>(
"lsst.ip.diffim.MaskedKernelSolution.build",
980 "Upper good pixel region: %d,%d -> %d,%d",
982 pexLog::TTrace<5>(
"lsst.ip.diffim.MaskedKernelSolution.build",
983 "Bottom good pixel region: %d,%d -> %d,%d",
985 pexLog::TTrace<5>(
"lsst.ip.diffim.MaskedKernelSolution.build",
986 "Left good pixel region: %d,%d -> %d,%d",
988 pexLog::TTrace<5>(
"lsst.ip.diffim.MaskedKernelSolution.build",
989 "Right good pixel region: %d,%d -> %d,%d",
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);
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();
1014 typename std::vector<afwGeom::Box2I>::iterator biter = boxArray.begin();
1015 for (; biter != boxArray.end(); ++biter) {
1016 int area = (*biter).getWidth() * (*biter).getHeight();
1024 Eigen::MatrixXd eiVarEstimate =
imageToEigenMatrix(sVarEstimate).array().inverse().matrix();
1026 eTemplate.resize(area, 1);
1027 eScience.resize(area, 1);
1028 eiVarEstimate.resize(area, 1);
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);
1039 std::vector<boost::shared_ptr<Eigen::MatrixXd> > convolvedEigenList(nKernelParameters);
1040 typename std::vector<boost::shared_ptr<Eigen::MatrixXd> >::iterator eiter =
1041 convolvedEigenList.
begin();
1043 for (kiter = basisList.begin(); kiter != basisList.end(); ++kiter, ++eiter) {
1045 Eigen::MatrixXd cMat(totalSize, 1);
1049 typename std::vector<afwGeom::Box2I>::iterator biter = boxArray.begin();
1050 for (; biter != boxArray.end(); ++biter) {
1051 int area = (*biter).getWidth() * (*biter).getHeight();
1055 esubimage.resize(area, 1);
1056 cMat.block(nTerms, 0, area, 1) = esubimage.block(0, 0, area, 1);
1061 boost::shared_ptr<Eigen::MatrixXd> cMatPtr (
new Eigen::MatrixXd(cMat));
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);
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);
1084 if (this->_fitForBackground)
1085 cMat.col(nParameters-1).fill(1.);
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)));
1093 new Eigen::MatrixXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_cMat))
1096 new Eigen::VectorXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_iVec))
1103 template <
typename InputT>
1106 bool fitForBackground,
1107 boost::shared_ptr<Eigen::MatrixXd> hMat,
1116 template <
typename InputT>
1118 Eigen::MatrixXd vMat = (this->_cMat)->jacobiSvd().matrixV();
1119 Eigen::MatrixXd vMatvMatT = vMat * vMat.transpose();
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) {
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);
1137 eValues(i) = 1.0 / eValues(i);
1140 Eigen::MatrixXd mInv = rMat * eValues.asDiagonal() * rMat.transpose();
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);
1154 Eigen::VectorXd term1 = (this->_aVec->transpose() * vMatvMatT * *(this->_aVec));
1155 if (term1.size() != 1)
1158 double term2a = (vMatvMatT * mLambda.inverse()).
trace();
1160 Eigen::VectorXd term2b = (this->_aVec->transpose() * (mInv * *(this->_bVec)));
1161 if (term2b.size() != 1)
1164 double risk = term1(0) + 2 * (term2a - term2b(0));
1165 pexLog::TTrace<6>(
"lsst.ip.diffim.RegularizedKernelSolution.estimateRisk",
1166 "Lambda = %.3f, Risk = %.5e",
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);
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]);
1178 return lambdas[index];
1182 template <
typename InputT>
1184 if (includeHmat ==
true) {
1185 return (boost::shared_ptr<Eigen::MatrixXd>(
1186 new Eigen::MatrixXd(*(this->_mMat) + _lambda * (*_hMat))
1194 template <
typename InputT>
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());
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;
1216 new Eigen::MatrixXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_cMat))
1219 new Eigen::VectorXd(this->_cMat->transpose() * this->_ivVec->asDiagonal() * *(this->_iVec))
1272 std::string lambdaType = _policy.getString(
"lambdaType");
1273 if (lambdaType ==
"absolute") {
1274 _lambda = _policy.getDouble(
"lambdaValue");
1276 else if (lambdaType ==
"relative") {
1277 _lambda = this->_mMat->trace() / this->_hMat->trace();
1278 _lambda *= _policy.getDouble(
"lambdaScaling");
1280 else if (lambdaType ==
"minimizeBiasedRisk") {
1281 double tol = _policy.getDouble(
"maxConditionNumber");
1282 _lambda = estimateRisk(tol);
1284 else if (lambdaType ==
"minimizeUnbiasedRisk") {
1285 _lambda = estimateRisk(std::numeric_limits<double>::max());
1291 pexLog::TTrace<5>(
"lsst.ip.diffim.KernelCandidate.build",
1292 "Applying kernel regularization with lambda = %.2e", _lambda);
1305 template <
typename InputT>
1307 std::vector<double> lambdas;
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);
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));
1341 _spatialKernelFunction(spatialKernelFunction),
1342 _constantFirstTerm(false),
1344 _background(background),
1352 bool isAlardLupton =
_policy.
getString(
"kernelBasisSet") ==
"alard-lupton";
1354 if (isAlardLupton || usePca) {
1369 boost::shared_ptr<Eigen::MatrixXd> mMat (
new Eigen::MatrixXd(
_nt,
_nt));
1370 boost::shared_ptr<Eigen::VectorXd> bVec (
new Eigen::VectorXd(
_nt));
1381 pexLog::TTrace<5>(
"lsst.ip.diffim.SpatialKernelSolution",
1382 "Initializing with size %d %d %d and constant first term = %s",
1389 boost::shared_ptr<Eigen::MatrixXd> qMat,
1390 boost::shared_ptr<Eigen::VectorXd> wVec) {
1392 pexLog::TTrace<8>(
"lsst.ip.diffim.SpatialKernelSolution.addConstraint",
1393 "Adding candidate at %f, %f", xCenter, yCenter);
1397 Eigen::VectorXd pK(
_nkt);
1399 for (
int idx = 0; idx <
_nkt; idx++) { paramsK[idx] = 0.0; }
1400 for (
int idx = 0; idx <
_nkt; idx++) {
1403 pK(idx) = (*_spatialKernelFunction)(xCenter, yCenter);
1406 Eigen::MatrixXd pKpKt = (pK * pK.transpose());
1409 Eigen::MatrixXd pBpBt;
1410 Eigen::MatrixXd pKpBt;
1412 pB = Eigen::VectorXd(
_nbt);
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++) {
1420 pB(idx) = (*_background)(xCenter, yCenter);
1423 pBpBt = (pB * pB.transpose());
1426 pKpBt = (pK * pB.transpose());
1430 std::cout <<
"Spatial weights" << std::endl;
1431 std::cout <<
"pKpKt " << pKpKt << std::endl;
1433 std::cout <<
"pBpBt " << pBpBt << std::endl;
1434 std::cout <<
"pKpBt " << pKpBt << std::endl;
1439 std::cout <<
"Spatial matrix inputs" << std::endl;
1440 std::cout <<
"M " << (*qMat) << std::endl;
1441 std::cout <<
"B " << (*wVec) << std::endl;
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();
1459 (*_bVec)(0) += (*wVec)(0);
1462 (*_mMat).block(0, mb, 1, _nbt) += (*qMat)(0,
_nbases) * pB.transpose();
1467 for(
int m1 = m0; m1 <
_nbases; m1++) {
1469 (*_mMat).block(m1*_nkt-dm, m1*_nkt-dm, _nkt, _nkt) +=
1470 (pKpKt * (*qMat)(m1,m1)).triangularView<Eigen::Upper>();
1473 for(
int m2 = m1+1; m2 <
_nbases; m2++) {
1474 (*_mMat).block(m1*_nkt-dm, m2*_nkt-dm, _nkt, _nkt) += (*qMat)(m1,m2) * pKpKt;
1479 (*_mMat).block(m1*_nkt-dm, mb, _nkt, _nbt) += (*qMat)(m1,
_nbases) * pKpBt;
1483 (*_bVec).segment(m1*_nkt-dm, _nkt) += (*wVec)(m1) * pK;
1488 (*_mMat).block(mb, mb, _nbt, _nbt) +=
1489 (pBpBt * (*qMat)(
_nbases,
_nbases)).triangularView<Eigen::Upper>();
1490 (*_bVec).segment(mb, _nbt) += (*wVec)(
_nbases) * pB;
1494 std::cout <<
"Spatial matrix outputs" << std::endl;
1495 std::cout <<
"mMat " << (*_mMat) << std::endl;
1496 std::cout <<
"bVec " << (*_bVec) << std::endl;
1508 (void)
_kernel->computeImage(*image,
false, pos[0], pos[1]);
1514 for (
int i = 0; i <
_nt; i++) {
1515 for (
int j = i+1; j <
_nt; j++) {
1516 (*_mMat)(j,i) = (*
_mMat)(i,j);
1547 std::vector<double> kCoeffs(
_nbases);
1548 for (
int i = 0; i <
_nbases; i++) {
1553 "I. Unable to determine spatial kernel solution %d (nan). Condition number = %.3e") % i % cNumber));
1555 kCoeffs[i] = (*_aVec)(i);
1560 new afwMath::LinearCombinationKernel(basisList, kCoeffs)
1566 std::vector<std::vector<double> > kCoeffs;
1568 for (
int i = 0, idx = 0; i <
_nbases; i++) {
1569 kCoeffs.push_back(std::vector<double>(
_nkt));
1577 "II. Unable to determine spatial kernel solution %d (nan). Condition number = %.3e") % idx % cNumber));
1579 kCoeffs[i][0] = (*_aVec)(idx++);
1582 for (
int j = 0; j <
_nkt; j++) {
1587 "III. Unable to determine spatial kernel solution %d (nan). Condition number = %.3e") % idx % cNumber));
1589 kCoeffs[i][j] = (*_aVec)(idx++);
1593 _kernel->setSpatialParameters(kCoeffs);
1603 for (
int i = 0; i <
_nbt; i++) {
1604 int idx =
_nt - _nbt + i;
1608 "Unable to determine spatial background solution %d (nan)") % (idx)));
1610 bgCoeffs[i] = (*_aVec)(idx);
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.
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
An include file to include the header files for lsst::afw::geom.
geom::Extent2I const getDimensions() const
Return the Kernel's dimensions (width, height)
Eigen matrix objects that present a view into an ndarray::Array.
Eigen::MatrixXi maskToEigenMatrix(lsst::afw::image::Mask< lsst::afw::image::MaskPixel > const &mask)
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.
Eigen3 view into an ndarray::Array.
virtual lsst::afw::image::Image< lsst::afw::math::Kernel::Pixel >::Ptr makeKernelImage()
bool getBool(const std::string &name) const
boost::uint16_t MaskPixel
a container for holding hierarchical configuration data in memory.
void setKernelParameters(std::vector< double > const ¶ms)
Set the kernel parameters of a spatially invariant kernel.
const std::string getString(const std::string &name) const
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
lsst::afw::math::Kernel::Ptr _kernel
Derived single-object convolution kernel.
Use (pixels & (given mask))
A Threshold is used to pass a threshold value to detection algorithms.
virtual lsst::afw::math::Kernel::Ptr getKernel()
geom::Box2I getBBox(ImageOrigin origin=PARENT) const
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
StaticKernelSolution(lsst::afw::math::KernelList const &basisList, bool fitForBackground)
boost::shared_ptr< Eigen::MatrixXd > getM()
virtual double getBackground()
An integer coordinate rectangle.
Declaration of classes to store the solution for convolution kernels.
table::Key< table::Array< Kernel::Pixel > > image
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)
An include file to include the header files for lsst::afw::image.
geom::Point2I getXY0() const
boost::shared_ptr< Kernel > Ptr
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's size; useful for passing to constructors.
An include file to include the header files for lsst::afw::detection.
lsst::afw::image::Image< lsst::afw::math::Kernel::Pixel > ImageT
afwMath::LinearCombinationKernel const & _kernel
virtual double getConditionNumber(ConditionNumberType conditionType)
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.
A kernel that is a linear combination of fixed basis kernels.
Vector< T, N > makeVector(T v1, T v2,..., T vN)
Variadic constructor for Vector.
boost::shared_ptr< Image< PixelT > > Ptr
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'd together.
int _nbt
Number of background terms.
bool _constantFirstTerm
Is the first term constant.
lsst::pex::policy::Policy _policy
Policy to control processing.
detail::SimpleInitializer< N > allocate(Vector< int, N > const &shape)
Create an expression that allocates uninitialized memory for an array.
#define LSST_EXCEPT(type,...)
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)
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'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.
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)
double estimateRisk(double maxCond)
#define LSST_EXCEPT_ADD(e, m)
A class to represent a 2-dimensional array of pixels.
std::vector< double > _createLambdaSteps()
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'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
int _nt
Total number of terms.