LSST Applications g0f08755f38+9c285cab97,g1635faa6d4+13f3999e92,g1653933729+a8ce1bb630,g1a0ca8cf93+bf6eb00ceb,g28da252d5a+0829b12dee,g29321ee8c0+5700dc9eac,g2bbee38e9b+9634bc57db,g2bc492864f+9634bc57db,g2cdde0e794+c2c89b37c4,g3156d2b45e+41e33cbcdc,g347aa1857d+9634bc57db,g35bb328faa+a8ce1bb630,g3a166c0a6a+9634bc57db,g3e281a1b8c+9f2c4e2fc3,g414038480c+077ccc18e7,g41af890bb2+fde0dd39b6,g5fbc88fb19+17cd334064,g781aacb6e4+a8ce1bb630,g80478fca09+55a9465950,g82479be7b0+d730eedb7d,g858d7b2824+9c285cab97,g9125e01d80+a8ce1bb630,g9726552aa6+10f999ec6a,ga5288a1d22+2a84bb7594,gacf8899fa4+c69c5206e8,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gc28159a63d+9634bc57db,gcf0d15dbbd+4b7d09cae4,gda3e153d99+9c285cab97,gda6a2b7d83+4b7d09cae4,gdaeeff99f8+1711a396fd,ge2409df99d+5e831397f4,ge79ae78c31+9634bc57db,gf0baf85859+147a0692ba,gf3967379c6+41c94011de,gf3fb38a9a8+8f07a9901b,gfb92a5be7c+9c285cab97,w.2024.46
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Protected Attributes | Static Protected Attributes | List of all members
lsst::ip::diffim::SpatialKernelSolution Class Reference

#include <KernelSolution.h>

Inheritance diagram for lsst::ip::diffim::SpatialKernelSolution:
lsst::ip::diffim::KernelSolution

Public Types

typedef std::shared_ptr< SpatialKernelSolutionPtr
 
enum  KernelSolvedBy {
  NONE = 0 , CHOLESKY_LDLT = 1 , CHOLESKY_LLT = 2 , LU = 3 ,
  EIGENVECTOR = 4
}
 
enum  ConditionNumberType { EIGENVALUE = 0 , SVD = 1 }
 
typedef lsst::afw::math::Kernel::Pixel PixelT
 
typedef lsst::afw::image::Image< lsst::afw::math::Kernel::PixelImageT
 

Public Member Functions

 SpatialKernelSolution (lsst::afw::math::KernelList const &basisList, lsst::afw::math::Kernel::SpatialFunctionPtr spatialKernelFunction, lsst::afw::math::Kernel::SpatialFunctionPtr background, lsst::daf::base::PropertySet const &ps)
 
virtual ~SpatialKernelSolution ()
 
void addConstraint (float xCenter, float yCenter, Eigen::MatrixXd const &qMat, Eigen::VectorXd const &wVec)
 
void solve ()
 
std::shared_ptr< lsst::afw::image::Image< lsst::afw::math::Kernel::Pixel > > makeKernelImage (lsst::geom::Point2D const &pos)
 
std::pair< std::shared_ptr< lsst::afw::math::LinearCombinationKernel >, lsst::afw::math::Kernel::SpatialFunctionPtrgetSolutionPair ()
 
virtual void solve (Eigen::MatrixXd const &mMat, Eigen::VectorXd const &bVec)
 
KernelSolvedBy getSolvedBy ()
 
virtual double getConditionNumber (ConditionNumberType conditionType)
 
virtual double getConditionNumber (Eigen::MatrixXd const &mMat, ConditionNumberType conditionType)
 
Eigen::MatrixXd const & getM ()
 
Eigen::VectorXd const & getB ()
 
void printM ()
 
void printB ()
 
void printA ()
 
int getId () const
 

Protected Attributes

int _id
 Unique ID for object.
 
Eigen::MatrixXd _mMat
 Derived least squares M matrix.
 
Eigen::VectorXd _bVec
 Derived least squares B vector.
 
Eigen::VectorXd _aVec
 Derived least squares solution matrix.
 
KernelSolvedBy _solvedBy
 Type of algorithm used to make solution.
 
bool _fitForBackground
 Background terms included in fit.
 

Static Protected Attributes

static int _SolutionId = 0
 Unique identifier for solution.
 

Detailed Description

Definition at line 174 of file KernelSolution.h.

Member Typedef Documentation

◆ ImageT

Definition at line 35 of file KernelSolution.h.

◆ PixelT

Definition at line 34 of file KernelSolution.h.

◆ Ptr

Definition at line 176 of file KernelSolution.h.

Member Enumeration Documentation

◆ ConditionNumberType

Enumerator
EIGENVALUE 
SVD 

Definition at line 45 of file KernelSolution.h.

◆ KernelSolvedBy

Enumerator
NONE 
CHOLESKY_LDLT 
CHOLESKY_LLT 
LU 
EIGENVECTOR 

Definition at line 37 of file KernelSolution.h.

Constructor & Destructor Documentation

◆ SpatialKernelSolution()

lsst::ip::diffim::SpatialKernelSolution::SpatialKernelSolution ( lsst::afw::math::KernelList const & basisList,
lsst::afw::math::Kernel::SpatialFunctionPtr spatialKernelFunction,
lsst::afw::math::Kernel::SpatialFunctionPtr background,
lsst::daf::base::PropertySet const & ps )

Definition at line 1283 of file KernelSolution.cc.

1288 :
1290 _spatialKernelFunction(spatialKernelFunction),
1291 _constantFirstTerm(false),
1292 _kernel(),
1293 _background(background),
1294 _kSum(0.0),
1295 _ps(ps.deepCopy()),
1296 _nbases(0),
1297 _nkt(0),
1298 _nbt(0),
1299 _nt(0) {
1300
1301 bool isAlardLupton = _ps->getAsString("kernelBasisSet") == "alard-lupton";
1302 bool usePca = _ps->getAsBool("usePcaForSpatialKernel");
1303 if (isAlardLupton || usePca) {
1304 _constantFirstTerm = true;
1305 }
1306 this->_fitForBackground = _ps->getAsBool("fitForBackground");
1307
1308 _nbases = basisList.size();
1309 _nkt = _spatialKernelFunction->getParameters().size();
1310 _nbt = _fitForBackground ? _background->getParameters().size() : 0;
1311 _nt = 0;
1312 if (_constantFirstTerm) {
1313 _nt = (_nbases - 1) * _nkt + 1 + _nbt;
1314 } else {
1315 _nt = _nbases * _nkt + _nbt;
1316 }
1317
1318 Eigen::MatrixXd mMat(_nt, _nt);
1319 Eigen::VectorXd bVec(_nt);
1320 mMat.setZero();
1321 bVec.setZero();
1322
1323 this->_mMat = mMat;
1324 this->_bVec = bVec;
1325
1327 new afwMath::LinearCombinationKernel(basisList, *_spatialKernelFunction)
1328 );
1329
1330 LOGL_DEBUG("TRACE3.ip.diffim.SpatialKernelSolution",
1331 "Initializing with size %d %d %d and constant first term = %s",
1332 _nkt, _nbt, _nt,
1333 _constantFirstTerm ? "true" : "false");
1334
1335 }
#define LOGL_DEBUG(logger, message...)
Log a debug-level message using a varargs/printf style interface.
Definition Log.h:515
std::vector< double > const & getParameters() const noexcept
Return all function parameters.
Definition Function.h:129
A kernel that is a linear combination of fixed basis kernels.
Definition Kernel.h:704
bool _fitForBackground
Background terms included in fit.
Eigen::VectorXd _bVec
Derived least squares B vector.
Eigen::MatrixXd _mMat
Derived least squares M matrix.
T size(T... args)

◆ ~SpatialKernelSolution()

virtual lsst::ip::diffim::SpatialKernelSolution::~SpatialKernelSolution ( )
inlinevirtual

Definition at line 185 of file KernelSolution.h.

185{};

Member Function Documentation

◆ addConstraint()

void lsst::ip::diffim::SpatialKernelSolution::addConstraint ( float xCenter,
float yCenter,
Eigen::MatrixXd const & qMat,
Eigen::VectorXd const & wVec )

Definition at line 1337 of file KernelSolution.cc.

1339 {
1340
1341 LOGL_DEBUG("TRACE5.ip.diffim.SpatialKernelSolution.addConstraint",
1342 "Adding candidate at %f, %f", xCenter, yCenter);
1343
1344 /* Calculate P matrices */
1345 /* Pure kernel terms */
1346 Eigen::VectorXd pK(_nkt);
1347 std::vector<double> paramsK = _spatialKernelFunction->getParameters();
1348 for (int idx = 0; idx < _nkt; idx++) { paramsK[idx] = 0.0; }
1349 for (int idx = 0; idx < _nkt; idx++) {
1350 paramsK[idx] = 1.0;
1351 _spatialKernelFunction->setParameters(paramsK);
1352 pK(idx) = (*_spatialKernelFunction)(xCenter, yCenter); /* Assume things don't vary over stamp */
1353 paramsK[idx] = 0.0;
1354 }
1355 Eigen::MatrixXd pKpKt = (pK * pK.transpose());
1356
1357 Eigen::VectorXd pB;
1358 Eigen::MatrixXd pBpBt;
1359 Eigen::MatrixXd pKpBt;
1360 if (_fitForBackground) {
1361 pB = Eigen::VectorXd(_nbt);
1362
1363 /* Pure background terms */
1364 std::vector<double> paramsB = _background->getParameters();
1365 for (int idx = 0; idx < _nbt; idx++) { paramsB[idx] = 0.0; }
1366 for (int idx = 0; idx < _nbt; idx++) {
1367 paramsB[idx] = 1.0;
1368 _background->setParameters(paramsB);
1369 pB(idx) = (*_background)(xCenter, yCenter); /* Assume things don't vary over stamp */
1370 paramsB[idx] = 0.0;
1371 }
1372 pBpBt = (pB * pB.transpose());
1373
1374 /* Cross terms */
1375 pKpBt = (pK * pB.transpose());
1376 }
1377
1378 if (DEBUG_MATRIX) {
1379 std::cout << "Spatial weights" << std::endl;
1380 std::cout << "pKpKt " << pKpKt << std::endl;
1381 if (_fitForBackground) {
1382 std::cout << "pBpBt " << pBpBt << std::endl;
1383 std::cout << "pKpBt " << pKpBt << std::endl;
1384 }
1385 }
1386
1387 if (DEBUG_MATRIX) {
1388 std::cout << "Spatial matrix inputs" << std::endl;
1389 std::cout << "M " << qMat << std::endl;
1390 std::cout << "B " << wVec << std::endl;
1391 }
1392
1393 /* first index to start the spatial blocks; default=0 for non-constant first term */
1394 int m0 = 0;
1395 /* how many rows/cols to adjust the matrices/vectors; default=0 for non-constant first term */
1396 int dm = 0;
1397 /* where to start the background terms; this is always true */
1398 int mb = _nt - _nbt;
1399
1400 if (_constantFirstTerm) {
1401 m0 = 1; /* we need to manually fill in the first (non-spatial) terms below */
1402 dm = _nkt-1; /* need to shift terms due to lack of spatial variation in first term */
1403
1404 _mMat(0, 0) += qMat(0,0);
1405 for(int m2 = 1; m2 < _nbases; m2++) {
1406 _mMat.block(0, m2*_nkt-dm, 1, _nkt) += qMat(0,m2) * pK.transpose();
1407 }
1408 _bVec(0) += wVec(0);
1409
1410 if (_fitForBackground) {
1411 _mMat.block(0, mb, 1, _nbt) += qMat(0,_nbases) * pB.transpose();
1412 }
1413 }
1414
1415 /* Fill in the spatial blocks */
1416 for(int m1 = m0; m1 < _nbases; m1++) {
1417 /* Diagonal kernel-kernel term; only use upper triangular part of pKpKt */
1418 _mMat.block(m1*_nkt-dm, m1*_nkt-dm, _nkt, _nkt) +=
1419 (pKpKt * qMat(m1,m1)).triangularView<Eigen::Upper>();
1420
1421 /* Kernel-kernel terms */
1422 for(int m2 = m1+1; m2 < _nbases; m2++) {
1423 _mMat.block(m1*_nkt-dm, m2*_nkt-dm, _nkt, _nkt) += qMat(m1,m2) * pKpKt;
1424 }
1425
1426 if (_fitForBackground) {
1427 /* Kernel cross terms with background */
1428 _mMat.block(m1*_nkt-dm, mb, _nkt, _nbt) += qMat(m1,_nbases) * pKpBt;
1429 }
1430
1431 /* B vector */
1432 _bVec.segment(m1*_nkt-dm, _nkt) += wVec(m1) * pK;
1433 }
1434
1435 if (_fitForBackground) {
1436 /* Background-background terms only */
1437 _mMat.block(mb, mb, _nbt, _nbt) +=
1438 (pBpBt * qMat(_nbases,_nbases)).triangularView<Eigen::Upper>();
1439 _bVec.segment(mb, _nbt) += wVec(_nbases) * pB;
1440 }
1441
1442 if (DEBUG_MATRIX) {
1443 std::cout << "Spatial matrix outputs" << std::endl;
1444 std::cout << "mMat " << _mMat << std::endl;
1445 std::cout << "bVec " << _bVec << std::endl;
1446 }
1447
1448 }
void setParameters(std::vector< double > const &params)
Set all function parameters.
Definition Function.h:156
T endl(T... args)
#define DEBUG_MATRIX

◆ getB()

Eigen::VectorXd const & lsst::ip::diffim::KernelSolution::getB ( )
inlineinherited

Definition at line 65 of file KernelSolution.h.

65{return _bVec;}

◆ getConditionNumber() [1/2]

double lsst::ip::diffim::KernelSolution::getConditionNumber ( ConditionNumberType conditionType)
virtualinherited

Definition at line 94 of file KernelSolution.cc.

94 {
95 return getConditionNumber(_mMat, conditionType);
96 }
virtual double getConditionNumber(ConditionNumberType conditionType)

◆ getConditionNumber() [2/2]

double lsst::ip::diffim::KernelSolution::getConditionNumber ( Eigen::MatrixXd const & mMat,
ConditionNumberType conditionType )
virtualinherited

Definition at line 98 of file KernelSolution.cc.

99 {
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 LOGL_DEBUG("TRACE3.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 LOGL_DEBUG("TRACE3.ip.diffim.KernelSolution.getConditionNumber",
118 "SVD eMax / eMin = %.3e", sMax / sMin);
119 return (sMax / sMin);
120 break;
121 }
122 default:
123 {
125 "Undefined ConditionNumberType : only EIGENVALUE, SVD allowed.");
126 break;
127 }
128 }
129 }
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition Exception.h:48
Reports invalid arguments.
Definition Runtime.h:66

◆ getId()

int lsst::ip::diffim::KernelSolution::getId ( ) const
inlineinherited

Definition at line 69 of file KernelSolution.h.

69{ return _id; }
int _id
Unique ID for object.

◆ getM()

Eigen::MatrixXd const & lsst::ip::diffim::KernelSolution::getM ( )
inlineinherited

Definition at line 64 of file KernelSolution.h.

64{return _mMat;}

◆ getSolutionPair()

std::pair< std::shared_ptr< afwMath::LinearCombinationKernel >, afwMath::Kernel::SpatialFunctionPtr > lsst::ip::diffim::SpatialKernelSolution::getSolutionPair ( )

Definition at line 1480 of file KernelSolution.cc.

1480 {
1482 throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return solution");
1483 }
1484
1485 return std::make_pair(_kernel, _background);
1486 }
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
Provides consistent interface for LSST exceptions.
Definition Exception.h:107
T make_pair(T... args)

◆ getSolvedBy()

KernelSolvedBy lsst::ip::diffim::KernelSolution::getSolvedBy ( )
inlineinherited

Definition at line 60 of file KernelSolution.h.

60{return _solvedBy;}

◆ makeKernelImage()

std::shared_ptr< lsst::afw::image::Image< lsst::afw::math::Kernel::Pixel > > lsst::ip::diffim::SpatialKernelSolution::makeKernelImage ( lsst::geom::Point2D const & pos)

Definition at line 1450 of file KernelSolution.cc.

1450 {
1452 throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return image");
1453 }
1456 );
1457 (void)_kernel->computeImage(*image, false, pos[0], pos[1]);
1458 return image;
1459 }
afw::table::Key< afw::table::Array< ImagePixelT > > image
A class to represent a 2-dimensional array of pixels.
Definition Image.h:51
lsst::geom::Extent2I const getDimensions() const
Return the Kernel's dimensions (width, height)
Definition Kernel.h:212
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:76

◆ printA()

void lsst::ip::diffim::KernelSolution::printA ( )
inlineinherited

Definition at line 68 of file KernelSolution.h.

Eigen::VectorXd _aVec
Derived least squares solution matrix.

◆ printB()

void lsst::ip::diffim::KernelSolution::printB ( )
inlineinherited

Definition at line 67 of file KernelSolution.h.

◆ printM()

void lsst::ip::diffim::KernelSolution::printM ( )
inlineinherited

Definition at line 66 of file KernelSolution.h.

◆ solve() [1/2]

void lsst::ip::diffim::SpatialKernelSolution::solve ( )
virtual

Reimplemented from lsst::ip::diffim::KernelSolution.

Definition at line 1461 of file KernelSolution.cc.

1461 {
1462 /* Fill in the other half of mMat */
1463 for (int i = 0; i < _nt; i++) {
1464 for (int j = i+1; j < _nt; j++) {
1465 _mMat(j,i) = _mMat(i,j);
1466 }
1467 }
1468
1469 try {
1471 } catch (pexExcept::Exception &e) {
1472 LSST_EXCEPT_ADD(e, "Unable to solve spatial kernel matrix");
1473 throw e;
1474 }
1475 /* Turn matrices into _kernel and _background */
1476 _setKernel();
1477 }
#define LSST_EXCEPT_ADD(e, m)
Add the current location and a message to an existing exception before rethrowing it.
Definition Exception.h:54

◆ solve() [2/2]

void lsst::ip::diffim::KernelSolution::solve ( Eigen::MatrixXd const & mMat,
Eigen::VectorXd const & bVec )
virtualinherited

Definition at line 131 of file KernelSolution.cc.

132 {
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::cpu_timer t;
144
145 LOGL_DEBUG("TRACE2.ip.diffim.KernelSolution.solve",
146 "Solving for kernel");
147 _solvedBy = LU;
148 Eigen::FullPivLU<Eigen::MatrixXd> lu(mMat);
149 if (lu.isInvertible()) {
150 aVec = lu.solve(bVec);
151 } else {
152 LOGL_DEBUG("TRACE3.ip.diffim.KernelSolution.solve",
153 "Unable to determine kernel via LU");
154 /* LAST RESORT */
155 try {
156
158 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(mMat);
159 Eigen::MatrixXd const& rMat = eVecValues.eigenvectors();
160 Eigen::VectorXd eValues = eVecValues.eigenvalues();
161
162 for (int i = 0; i != eValues.rows(); ++i) {
163 if (eValues(i) != 0.0) {
164 eValues(i) = 1.0/eValues(i);
165 }
166 }
167
168 aVec = rMat * eValues.asDiagonal() * rMat.transpose() * bVec;
169 } catch (pexExcept::Exception& e) {
170
171 _solvedBy = NONE;
172 LOGL_DEBUG("TRACE3.ip.diffim.KernelSolution.solve",
173 "Unable to determine kernel via eigen-values");
174
175 throw LSST_EXCEPT(pexExcept::Exception, "Unable to determine kernel solution");
176 }
177 }
178
179 t.stop();
180 double time = 1e-9 * t.elapsed().wall;
181 LOGL_DEBUG("TRACE3.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 = aVec;
190 }

Member Data Documentation

◆ _aVec

Eigen::VectorXd lsst::ip::diffim::KernelSolution::_aVec
protectedinherited

Derived least squares solution matrix.

Definition at line 75 of file KernelSolution.h.

◆ _bVec

Eigen::VectorXd lsst::ip::diffim::KernelSolution::_bVec
protectedinherited

Derived least squares B vector.

Definition at line 74 of file KernelSolution.h.

◆ _fitForBackground

bool lsst::ip::diffim::KernelSolution::_fitForBackground
protectedinherited

Background terms included in fit.

Definition at line 77 of file KernelSolution.h.

◆ _id

int lsst::ip::diffim::KernelSolution::_id
protectedinherited

Unique ID for object.

Definition at line 72 of file KernelSolution.h.

◆ _mMat

Eigen::MatrixXd lsst::ip::diffim::KernelSolution::_mMat
protectedinherited

Derived least squares M matrix.

Definition at line 73 of file KernelSolution.h.

◆ _SolutionId

int lsst::ip::diffim::KernelSolution::_SolutionId = 0
staticprotectedinherited

Unique identifier for solution.

Definition at line 78 of file KernelSolution.h.

◆ _solvedBy

KernelSolvedBy lsst::ip::diffim::KernelSolution::_solvedBy
protectedinherited

Type of algorithm used to make solution.

Definition at line 76 of file KernelSolution.h.


The documentation for this class was generated from the following files: