LSSTApplications  18.1.0
LSSTDataManagementBasePackage
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::pex::policy::Policy policy)
 
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::afw::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. More...
 
Eigen::MatrixXd _mMat
 Derived least squares M matrix. More...
 
Eigen::VectorXd _bVec
 Derived least squares B vector. More...
 
Eigen::VectorXd _aVec
 Derived least squares solution matrix. More...
 
KernelSolvedBy _solvedBy
 Type of algorithm used to make solution. More...
 
bool _fitForBackground
 Background terms included in fit. More...
 

Static Protected Attributes

static int _SolutionId = 0
 Unique identifier for solution. More...
 

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

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::pex::policy::Policy  policy 
)

Definition at line 1281 of file KernelSolution.cc.

1286  :
1287  KernelSolution(),
1288  _spatialKernelFunction(spatialKernelFunction),
1289  _constantFirstTerm(false),
1290  _kernel(),
1291  _background(background),
1292  _kSum(0.0),
1293  _policy(policy),
1294  _nbases(0),
1295  _nkt(0),
1296  _nbt(0),
1297  _nt(0) {
1298 
1299  bool isAlardLupton = _policy.getString("kernelBasisSet") == "alard-lupton";
1300  bool usePca = _policy.getBool("usePcaForSpatialKernel");
1301  if (isAlardLupton || usePca) {
1302  _constantFirstTerm = true;
1303  }
1304  this->_fitForBackground = _policy.getBool("fitForBackground");
1305 
1306  _nbases = basisList.size();
1307  _nkt = _spatialKernelFunction->getParameters().size();
1308  _nbt = _fitForBackground ? _background->getParameters().size() : 0;
1309  _nt = 0;
1310  if (_constantFirstTerm) {
1311  _nt = (_nbases - 1) * _nkt + 1 + _nbt;
1312  } else {
1313  _nt = _nbases * _nkt + _nbt;
1314  }
1315 
1316  Eigen::MatrixXd mMat(_nt, _nt);
1317  Eigen::VectorXd bVec(_nt);
1318  mMat.setZero();
1319  bVec.setZero();
1320 
1321  this->_mMat = mMat;
1322  this->_bVec = bVec;
1323 
1325  new afwMath::LinearCombinationKernel(basisList, *_spatialKernelFunction)
1326  );
1327 
1328  LOGL_DEBUG("TRACE3.ip.diffim.SpatialKernelSolution",
1329  "Initializing with size %d %d %d and constant first term = %s",
1330  _nkt, _nbt, _nt,
1331  _constantFirstTerm ? "true" : "false");
1332 
1333  }
Eigen::VectorXd _bVec
Derived least squares B vector.
const std::string getString(const std::string &name) const
return a string value associated with the given name .
Definition: Policy.h:631
#define LOGL_DEBUG(logger, message...)
Log a debug-level message using a varargs/printf style interface.
Definition: Log.h:489
bool _fitForBackground
Background terms included in fit.
A kernel that is a linear combination of fixed basis kernels.
Definition: Kernel.h:741
bool getBool(const std::string &name) const
return a boolean value associated with the given name.
Definition: Policy.h:589
Eigen::MatrixXd _mMat
Derived least squares M matrix.

◆ ~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 1335 of file KernelSolution.cc.

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

◆ getB()

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

Definition at line 65 of file KernelSolution.h.

65 {return _bVec;}
Eigen::VectorXd _bVec
Derived least squares B vector.

◆ getConditionNumber() [1/2]

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

Definition at line 92 of file KernelSolution.cc.

92  {
93  return getConditionNumber(_mMat, conditionType);
94  }
virtual double getConditionNumber(ConditionNumberType conditionType)
Eigen::MatrixXd _mMat
Derived least squares M matrix.

◆ getConditionNumber() [2/2]

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

Definition at line 96 of file KernelSolution.cc.

97  {
98  switch (conditionType) {
99  case EIGENVALUE:
100  {
101  Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(mMat);
102  Eigen::VectorXd eValues = eVecValues.eigenvalues();
103  double eMax = eValues.maxCoeff();
104  double eMin = eValues.minCoeff();
105  LOGL_DEBUG("TRACE3.ip.diffim.KernelSolution.getConditionNumber",
106  "EIGENVALUE eMax / eMin = %.3e", eMax / eMin);
107  return (eMax / eMin);
108  break;
109  }
110  case SVD:
111  {
112  Eigen::VectorXd sValues = mMat.jacobiSvd().singularValues();
113  double sMax = sValues.maxCoeff();
114  double sMin = sValues.minCoeff();
115  LOGL_DEBUG("TRACE3.ip.diffim.KernelSolution.getConditionNumber",
116  "SVD eMax / eMin = %.3e", sMax / sMin);
117  return (sMax / sMin);
118  break;
119  }
120  default:
121  {
123  "Undefined ConditionNumberType : only EIGENVALUE, SVD allowed.");
124  break;
125  }
126  }
127  }
#define LOGL_DEBUG(logger, message...)
Log a debug-level message using a varargs/printf style interface.
Definition: Log.h:489
#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;}
Eigen::MatrixXd _mMat
Derived least squares M matrix.

◆ getSolutionPair()

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

Definition at line 1478 of file KernelSolution.cc.

1478  {
1480  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return solution");
1481  }
1482 
1483  return std::make_pair(_kernel, _background);
1484  }
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
Provides consistent interface for LSST exceptions.
Definition: Exception.h:107
T make_pair(T... args)
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48

◆ getSolvedBy()

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

Definition at line 60 of file KernelSolution.h.

60 {return _solvedBy;}
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.

◆ makeKernelImage()

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

Definition at line 1448 of file KernelSolution.cc.

1448  {
1450  throw LSST_EXCEPT(pexExcept::Exception, "Kernel not solved; cannot return image");
1451  }
1453  new afwImage::Image<afwMath::Kernel::Pixel>(_kernel->getDimensions())
1454  );
1455  (void)_kernel->computeImage(*image, false, pos[0], pos[1]);
1456  return image;
1457  }
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
Provides consistent interface for LSST exceptions.
Definition: Exception.h:107
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
afw::table::Key< afw::table::Array< ImagePixelT > > image
Backwards-compatibility support for depersisting the old Calib (FluxMag0/FluxMag0Err) objects...
A class to represent a 2-dimensional array of pixels.
Definition: Image.h:59

◆ printA()

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

Definition at line 68 of file KernelSolution.h.

68 {std::cout << _aVec << std::endl;}
T endl(T... args)
Eigen::VectorXd _aVec
Derived least squares solution matrix.

◆ printB()

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

Definition at line 67 of file KernelSolution.h.

67 {std::cout << _bVec << std::endl;}
Eigen::VectorXd _bVec
Derived least squares B vector.
T endl(T... args)

◆ printM()

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

Definition at line 66 of file KernelSolution.h.

66 {std::cout << _mMat << std::endl;}
T endl(T... args)
Eigen::MatrixXd _mMat
Derived least squares M matrix.

◆ solve() [1/2]

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

Definition at line 129 of file KernelSolution.cc.

130  {
131 
132  if (DEBUG_MATRIX) {
133  std::cout << "M " << std::endl;
134  std::cout << mMat << std::endl;
135  std::cout << "B " << std::endl;
136  std::cout << bVec << std::endl;
137  }
138 
139  Eigen::VectorXd aVec = Eigen::VectorXd::Zero(bVec.size());
140 
141  boost::timer t;
142  t.restart();
143 
144  LOGL_DEBUG("TRACE2.ip.diffim.KernelSolution.solve",
145  "Solving for kernel");
146  _solvedBy = LU;
147  Eigen::FullPivLU<Eigen::MatrixXd> lu(mMat);
148  if (lu.isInvertible()) {
149  aVec = lu.solve(bVec);
150  } else {
151  LOGL_DEBUG("TRACE3.ip.diffim.KernelSolution.solve",
152  "Unable to determine kernel via LU");
153  /* LAST RESORT */
154  try {
155 
157  Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(mMat);
158  Eigen::MatrixXd const& rMat = eVecValues.eigenvectors();
159  Eigen::VectorXd eValues = eVecValues.eigenvalues();
160 
161  for (int i = 0; i != eValues.rows(); ++i) {
162  if (eValues(i) != 0.0) {
163  eValues(i) = 1.0/eValues(i);
164  }
165  }
166 
167  aVec = rMat * eValues.asDiagonal() * rMat.transpose() * bVec;
168  } catch (pexExcept::Exception& e) {
169 
170  _solvedBy = NONE;
171  LOGL_DEBUG("TRACE3.ip.diffim.KernelSolution.solve",
172  "Unable to determine kernel via eigen-values");
173 
174  throw LSST_EXCEPT(pexExcept::Exception, "Unable to determine kernel solution");
175  }
176  }
177 
178  double time = t.elapsed();
179  LOGL_DEBUG("TRACE3.ip.diffim.KernelSolution.solve",
180  "Compute time for matrix math : %.2f s", time);
181 
182  if (DEBUG_MATRIX) {
183  std::cout << "A " << std::endl;
184  std::cout << aVec << std::endl;
185  }
186 
187  _aVec = aVec;
188  }
T endl(T... args)
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
Provides consistent interface for LSST exceptions.
Definition: Exception.h:107
#define LOGL_DEBUG(logger, message...)
Log a debug-level message using a varargs/printf style interface.
Definition: Log.h:489
T time(T... args)
#define DEBUG_MATRIX
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
Eigen::VectorXd _aVec
Derived least squares solution matrix.

◆ solve() [2/2]

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

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

Definition at line 1459 of file KernelSolution.cc.

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

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: