#include <KernelSolution.h>
Definition at line 31 of file KernelSolution.h.
◆ ImageT
◆ PixelT
◆ Ptr
◆ ConditionNumberType
◆ KernelSolvedBy
Enumerator |
---|
NONE | |
CHOLESKY_LDLT | |
CHOLESKY_LLT | |
LU | |
EIGENVECTOR | |
Definition at line 37 of file KernelSolution.h.
◆ KernelSolution() [1/3]
lsst::ip::diffim::KernelSolution::KernelSolution |
( |
Eigen::MatrixXd |
mMat, |
|
|
Eigen::VectorXd |
bVec, |
|
|
bool |
fitForBackground |
|
) |
| |
|
explicit |
Definition at line 55 of file KernelSolution.cc.
Eigen::VectorXd _bVec
Derived least squares B vector.
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
static int _SolutionId
Unique identifier for solution.
bool _fitForBackground
Background terms included in fit.
int _id
Unique ID for object.
Eigen::MatrixXd _mMat
Derived least squares M matrix.
Eigen::VectorXd _aVec
Derived least squares solution matrix.
◆ KernelSolution() [2/3]
lsst::ip::diffim::KernelSolution::KernelSolution |
( |
bool |
fitForBackground | ) |
|
|
explicit |
Definition at line 68 of file KernelSolution.cc.
Eigen::VectorXd _bVec
Derived least squares B vector.
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
static int _SolutionId
Unique identifier for solution.
bool _fitForBackground
Background terms included in fit.
int _id
Unique ID for object.
Eigen::MatrixXd _mMat
Derived least squares M matrix.
Eigen::VectorXd _aVec
Derived least squares solution matrix.
◆ KernelSolution() [3/3]
lsst::ip::diffim::KernelSolution::KernelSolution |
( |
| ) |
|
|
explicit |
Definition at line 79 of file KernelSolution.cc.
Eigen::VectorXd _bVec
Derived least squares B vector.
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
static int _SolutionId
Unique identifier for solution.
bool _fitForBackground
Background terms included in fit.
int _id
Unique ID for object.
Eigen::MatrixXd _mMat
Derived least squares M matrix.
Eigen::VectorXd _aVec
Derived least squares solution matrix.
◆ ~KernelSolution()
virtual lsst::ip::diffim::KernelSolution::~KernelSolution |
( |
| ) |
|
|
inlinevirtual |
◆ getB()
Eigen::VectorXd const& lsst::ip::diffim::KernelSolution::getB |
( |
| ) |
|
|
inline |
Definition at line 65 of file KernelSolution.h.
Eigen::VectorXd _bVec
Derived least squares B vector.
◆ getConditionNumber() [1/2]
double lsst::ip::diffim::KernelSolution::getConditionNumber |
( |
ConditionNumberType |
conditionType | ) |
|
|
virtual |
Definition at line 92 of file KernelSolution.cc.
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 |
|
) |
| |
|
virtual |
Definition at line 96 of file KernelSolution.cc.
98 switch (conditionType) {
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);
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);
123 "Undefined ConditionNumberType : only EIGENVALUE, SVD allowed.");
#define LOGL_DEBUG(logger, message...)
Log a debug-level message using a varargs/printf style interface.
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Reports invalid arguments.
◆ getId()
int lsst::ip::diffim::KernelSolution::getId |
( |
| ) |
const |
|
inline |
◆ getM()
Eigen::MatrixXd const& lsst::ip::diffim::KernelSolution::getM |
( |
| ) |
|
|
inline |
Definition at line 64 of file KernelSolution.h.
Eigen::MatrixXd _mMat
Derived least squares M matrix.
◆ getSolvedBy()
Definition at line 60 of file KernelSolution.h.
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
◆ printA()
void lsst::ip::diffim::KernelSolution::printA |
( |
| ) |
|
|
inline |
Definition at line 68 of file KernelSolution.h.
Eigen::VectorXd _aVec
Derived least squares solution matrix.
◆ printB()
void lsst::ip::diffim::KernelSolution::printB |
( |
| ) |
|
|
inline |
Definition at line 67 of file KernelSolution.h.
Eigen::VectorXd _bVec
Derived least squares B vector.
◆ printM()
void lsst::ip::diffim::KernelSolution::printM |
( |
| ) |
|
|
inline |
Definition at line 66 of file KernelSolution.h.
Eigen::MatrixXd _mMat
Derived least squares M matrix.
◆ solve() [1/2]
void lsst::ip::diffim::KernelSolution::solve |
( |
| ) |
|
|
virtual |
◆ solve() [2/2]
void lsst::ip::diffim::KernelSolution::solve |
( |
Eigen::MatrixXd const & |
mMat, |
|
|
Eigen::VectorXd const & |
bVec |
|
) |
| |
|
virtual |
Definition at line 129 of file KernelSolution.cc.
139 Eigen::VectorXd aVec = Eigen::VectorXd::Zero(bVec.size());
144 LOGL_DEBUG(
"TRACE2.ip.diffim.KernelSolution.solve",
145 "Solving for kernel");
147 Eigen::FullPivLU<Eigen::MatrixXd> lu(mMat);
148 if (lu.isInvertible()) {
149 aVec = lu.solve(bVec);
151 LOGL_DEBUG(
"TRACE3.ip.diffim.KernelSolution.solve",
152 "Unable to determine kernel via LU");
157 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(mMat);
158 Eigen::MatrixXd
const& rMat = eVecValues.eigenvectors();
159 Eigen::VectorXd eValues = eVecValues.eigenvalues();
161 for (
int i = 0; i != eValues.rows(); ++i) {
162 if (eValues(i) != 0.0) {
163 eValues(i) = 1.0/eValues(i);
167 aVec = rMat * eValues.asDiagonal() * rMat.transpose() * bVec;
171 LOGL_DEBUG(
"TRACE3.ip.diffim.KernelSolution.solve",
172 "Unable to determine kernel via eigen-values");
178 double time = t.elapsed();
179 LOGL_DEBUG(
"TRACE3.ip.diffim.KernelSolution.solve",
180 "Compute time for matrix math : %.2f s", time);
KernelSolvedBy _solvedBy
Type of algorithm used to make solution.
Provides consistent interface for LSST exceptions.
#define LOGL_DEBUG(logger, message...)
Log a debug-level message using a varargs/printf style interface.
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Eigen::VectorXd _aVec
Derived least squares solution matrix.
◆ _aVec
Eigen::VectorXd lsst::ip::diffim::KernelSolution::_aVec |
|
protected |
◆ _bVec
Eigen::VectorXd lsst::ip::diffim::KernelSolution::_bVec |
|
protected |
◆ _fitForBackground
bool lsst::ip::diffim::KernelSolution::_fitForBackground |
|
protected |
◆ _id
int lsst::ip::diffim::KernelSolution::_id |
|
protected |
◆ _mMat
Eigen::MatrixXd lsst::ip::diffim::KernelSolution::_mMat |
|
protected |
◆ _SolutionId
int lsst::ip::diffim::KernelSolution::_SolutionId = 0 |
|
staticprotected |
◆ _solvedBy
The documentation for this class was generated from the following files:
- /j/snowflake/release/lsstsw/stack/Linux64/ip_diffim/18.1.0/include/lsst/ip/diffim/KernelSolution.h
- /j/snowflake/release/lsstsw/stack/Linux64/ip_diffim/18.1.0/src/KernelSolution.cc