LSST Applications g0b6bd0c080+a72a5dd7e6,g1182afd7b4+2a019aa3bb,g17e5ecfddb+2b8207f7de,g1d67935e3f+06cf436103,g38293774b4+ac198e9f13,g396055baef+6a2097e274,g3b44f30a73+6611e0205b,g480783c3b1+98f8679e14,g48ccf36440+89c08d0516,g4b93dc025c+98f8679e14,g5c4744a4d9+a302e8c7f0,g613e996a0d+e1c447f2e0,g6c8d09e9e7+25247a063c,g7271f0639c+98f8679e14,g7a9cd813b8+124095ede6,g9d27549199+a302e8c7f0,ga1cf026fa3+ac198e9f13,ga32aa97882+7403ac30ac,ga786bb30fb+7a139211af,gaa63f70f4e+9994eb9896,gabf319e997+ade567573c,gba47b54d5d+94dc90c3ea,gbec6a3398f+06cf436103,gc6308e37c7+07dd123edb,gc655b1545f+ade567573c,gcc9029db3c+ab229f5caf,gd01420fc67+06cf436103,gd877ba84e5+06cf436103,gdb4cecd868+6f279b5b48,ge2d134c3d5+cc4dbb2e3f,ge448b5faa6+86d1ceac1d,gecc7e12556+98f8679e14,gf3ee170dca+25247a063c,gf4ac96e456+ade567573c,gf9f5ea5b4d+ac198e9f13,gff490e6085+8c2580be5c,w.2022.27
LSST Data Management Base Package
ImagePca.cc
Go to the documentation of this file.
1// -*- lsst-c++ -*-
2
3/*
4 * LSST Data Management System
5 * Copyright 2008, 2009, 2010 LSST Corporation.
6 *
7 * This product includes software developed by the
8 * LSST Project (http://www.lsst.org/).
9 *
10 * This program is free software: you can redistribute it and/or modify
11 * it under the terms of the GNU General Public License as published by
12 * the Free Software Foundation, either version 3 of the License, or
13 * (at your option) any later version.
14 *
15 * This program is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the LSST License Statement and
21 * the GNU General Public License along with this program. If not,
22 * see <http://www.lsstcorp.org/LegalNotices/>.
23 */
24
25/*
26 * Utilities to support PCA analysis of a set of images
27 */
28#include <algorithm>
29#include <cmath>
30#include <cstdint>
31#include <memory>
32
33#include "Eigen/Core"
34#include "Eigen/Eigenvalues"
35
38
39namespace afwMath = lsst::afw::math;
40
41namespace lsst {
42namespace afw {
43namespace image {
44
45template <typename ImageT>
46ImagePca<ImageT>::ImagePca(bool constantWeight)
47 : _imageList(),
48 _fluxList(),
49 _dimensions(0, 0),
50 _constantWeight(constantWeight),
51 _eigenValues(std::vector<double>()),
52 _eigenImages(ImageList()) {}
53
54template <typename ImageT>
55ImagePca<ImageT>::ImagePca(ImagePca const&) = default;
56template <typename ImageT>
58template <typename ImageT>
60template <typename ImageT>
62
63template <typename ImageT>
65 if (_imageList.empty()) {
66 _dimensions = img->getDimensions();
67 } else {
68 if (getDimensions() != img->getDimensions()) {
70 (boost::format("Dimension mismatch: saw %dx%d; expected %dx%d") %
71 img->getWidth() % img->getHeight() % _dimensions.getX() % _dimensions.getY())
72 .str());
73 }
74 }
75
76 if (flux == 0.0) {
77 throw LSST_EXCEPT(lsst::pex::exceptions::OutOfRangeError, "Flux may not be zero");
78 }
79
80 _imageList.push_back(img);
81 _fluxList.push_back(flux);
82}
83
84template <typename ImageT>
86 return _imageList;
87}
88
89template <typename ImageT>
91 if (_imageList.empty()) {
92 throw LSST_EXCEPT(lsst::pex::exceptions::LengthError, "You haven't provided any images");
93 }
94
95 std::shared_ptr<ImageT> mean(new ImageT(getDimensions()));
96 *mean = static_cast<typename ImageT::Pixel>(0);
97
98 for (typename ImageList::const_iterator ptr = _imageList.begin(), end = _imageList.end(); ptr != end;
99 ++ptr) {
100 *mean += **ptr;
101 }
102 *mean /= _imageList.size();
103
104 return mean;
105}
106
107namespace {
108/*
109 * Analyze the images in an ImagePca, calculating the PCA decomposition (== Karhunen-Lo\`eve basis)
110 *
111 * The notation is that in chapter 7 of Gyula Szokoly's thesis at JHU
112 */
113template <typename T>
114struct SortEvalueDecreasing
115 : public std::binary_function<std::pair<T, int> const&, std::pair<T, int> const&, bool> {
116 bool operator()(std::pair<T, int> const& a, std::pair<T, int> const& b) {
117 return a.first > b.first; // N.b. sort on greater
118 }
119};
120} // namespace
121
122template <typename ImageT>
124 int const nImage = _imageList.size();
125
126 if (nImage == 0) {
127 throw LSST_EXCEPT(lsst::pex::exceptions::LengthError, "No images provided for PCA analysis");
128 }
129 /*
130 * Eigen doesn't like 1x1 matrices, but we don't really need it to handle a single matrix...
131 */
132 if (nImage == 1) {
133 _eigenImages.clear();
134 _eigenImages.push_back(std::shared_ptr<ImageT>(new ImageT(*_imageList[0], true)));
135
136 _eigenValues.clear();
137 _eigenValues.push_back(1.0);
138
139 return;
140 }
141 /*
142 * Find the eigenvectors/values of the scalar product matrix, R' (Eq. 7.4)
143 */
144 Eigen::MatrixXd R(nImage, nImage); // residuals' inner products
145
146 double flux_bar = 0; // mean of flux for all regions
147 for (int i = 0; i != nImage; ++i) {
148 typename GetImage<ImageT>::type const& im_i = *GetImage<ImageT>::getImage(_imageList[i]);
149 double const flux_i = getFlux(i);
150 flux_bar += flux_i;
151
152 for (int j = i; j != nImage; ++j) {
153 typename GetImage<ImageT>::type const& im_j = *GetImage<ImageT>::getImage(_imageList[j]);
154 double const flux_j = getFlux(j);
155
156 double dot = innerProduct(im_i, im_j);
157 if (_constantWeight) {
158 dot /= flux_i * flux_j;
159 }
160 R(i, j) = R(j, i) = dot / nImage;
161 }
162 }
163 flux_bar /= nImage;
164 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(R);
165 Eigen::MatrixXd const& Q = eVecValues.eigenvectors();
166 Eigen::VectorXd const& lambda = eVecValues.eigenvalues();
167 //
168 // We need to sort the eigenValues, and remember the permutation we applied to the eigenImages
169 // We'll use the vector lambdaAndIndex to achieve this
170 //
171 std::vector<std::pair<double, int> > lambdaAndIndex; // pairs (eValue, index)
172 lambdaAndIndex.reserve(nImage);
173
174 for (int i = 0; i != nImage; ++i) {
175 lambdaAndIndex.emplace_back(lambda(i), i);
176 }
177 std::sort(lambdaAndIndex.begin(), lambdaAndIndex.end(), SortEvalueDecreasing<double>());
178 //
179 // Save the (sorted) eigen values
180 //
181 _eigenValues.clear();
182 _eigenValues.reserve(nImage);
183 for (int i = 0; i != nImage; ++i) {
184 _eigenValues.push_back(lambdaAndIndex[i].first);
185 }
186 //
187 // Contruct the first ncomp eigenimages in basis
188 //
189 int ncomp = 100; // number of components to keep
190
191 _eigenImages.clear();
192 _eigenImages.reserve(ncomp < nImage ? ncomp : nImage);
193
194 for (int i = 0; i < ncomp; ++i) {
195 if (i >= nImage) {
196 continue;
197 }
198
199 int const ii = lambdaAndIndex[i].second; // the index after sorting (backwards) by eigenvalue
200
201 std::shared_ptr<ImageT> eImage(new ImageT(_dimensions));
202 *eImage = static_cast<typename ImageT::Pixel>(0);
203
204 for (int j = 0; j != nImage; ++j) {
205 int const jj = lambdaAndIndex[j].second; // the index after sorting (backwards) by eigenvalue
206 double const weight = Q(jj, ii) * (_constantWeight ? flux_bar / getFlux(jj) : 1);
207 eImage->scaledPlus(weight, *_imageList[jj]);
208 }
209 _eigenImages.push_back(eImage);
210 }
211}
212
213namespace {
214/*
215 * Fit a LinearCombinationKernel to an Image, allowing the coefficients of the components to vary
216 *
217 * return std::pair(best-fit kernel, std::pair(amp, chi^2))
218 */
219template <typename MaskedImageT>
221 typename ImagePca<MaskedImageT>::ImageList const& eigenImages, // Eigen images
222 int nEigen, // Number of eigen images to use
223 MaskedImageT const& image // The image to be fit
224) {
225 using ImageT = typename MaskedImageT::Image;
226
227 if (nEigen == 0) {
228 throw LSST_EXCEPT(lsst::pex::exceptions::LengthError, "You must have at least one eigen image");
229 } else if (nEigen > static_cast<int>(eigenImages.size())) {
231 (boost::format("You only have %d eigen images (you asked for %d)") %
232 eigenImages.size() % nEigen)
233 .str());
234 }
235 /*
236 * Solve the linear problem image = sum x_i K_i + epsilon; we solve this for x_i by constructing the
237 * normal equations, A x = b
238 */
239 Eigen::MatrixXd A(nEigen, nEigen);
240 Eigen::VectorXd b(nEigen);
241
242 for (int i = 0; i != nEigen; ++i) {
243 b(i) = innerProduct(*eigenImages[i]->getImage(), *image.getImage());
244
245 for (int j = i; j != nEigen; ++j) {
246 A(i, j) = A(j, i) = innerProduct(*eigenImages[i]->getImage(), *eigenImages[j]->getImage());
247 }
248 }
249 Eigen::VectorXd x(nEigen);
250
251 if (nEigen == 1) {
252 x(0) = b(0) / A(0, 0);
253 } else {
254 x = A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
255 }
256 //
257 // Accumulate the best-fit-image in bestFitImage
258 //
259 std::shared_ptr<ImageT> bestFitImage = std::make_shared<ImageT>(eigenImages[0]->getDimensions());
260
261 for (int i = 0; i != nEigen; ++i) {
262 bestFitImage->scaledPlus(x[i], *eigenImages[i]->getImage());
263 }
264
265 return bestFitImage;
266}
267
268template <typename ImageT>
269double do_updateBadPixels(detail::basic_tag const&, typename ImagePca<ImageT>::ImageList const&,
270 std::vector<double> const&, typename ImagePca<ImageT>::ImageList const&,
271 unsigned long, int const) {
272 return 0.0;
273}
274
275template <typename ImageT>
276double do_updateBadPixels(detail::MaskedImage_tag const&,
277 typename ImagePca<ImageT>::ImageList const& imageList,
278 std::vector<double> const& fluxes, // fluxes of images
279 typename ImagePca<ImageT>::ImageList const& eigenImages, // Eigen images
280 unsigned long mask,
281 int const ncomp
282) {
283 int const nImage = imageList.size();
284 if (nImage == 0) {
286 "Please provide at least one Image for me to update");
287 }
288 lsst::geom::Extent2I dimensions = imageList[0]->getDimensions();
289 int const height = dimensions.getY();
290
291 double maxChange = 0.0; // maximum change to the input images
292
293 if (ncomp == 0) { // use mean of good pixels
294 typename ImageT::Image mean(dimensions); // desired mean image
295 image::Image<float> weight(mean.getDimensions()); // weight of each pixel
296
297 for (int i = 0; i != nImage; ++i) {
298 double const flux_i = fluxes[i];
299
300 for (int y = 0; y != height; ++y) {
301 typename ImageT::const_x_iterator iptr = imageList[i]->row_begin(y);
302 image::Image<float>::x_iterator wptr = weight.row_begin(y);
303 for (typename ImageT::Image::x_iterator mptr = mean.row_begin(y), end = mean.row_end(y);
304 mptr != end; ++mptr, ++iptr, ++wptr) {
305 if (!(iptr.mask() & mask) && iptr.variance() > 0.0) {
306 typename ImageT::Image::Pixel value = iptr.image() / flux_i;
307 float const var = iptr.variance() / (flux_i * flux_i);
308 float const ivar = 1.0 / var;
309
310 if (std::isfinite(value * ivar)) {
311 *mptr += value * ivar;
312 *wptr += ivar;
313 }
314 }
315 }
316 }
317 }
318 //
319 // Calculate mean
320 //
321 for (int y = 0; y != height; ++y) {
322 image::Image<float>::x_iterator wptr = weight.row_begin(y);
323 for (typename ImageT::Image::x_iterator mptr = mean.row_begin(y), end = mean.row_end(y);
324 mptr != end; ++mptr, ++wptr) {
325 if (*wptr == 0.0) { // oops, no good values
326 ;
327 } else {
328 *mptr /= *wptr;
329 }
330 }
331 }
332 //
333 // Replace bad values by mean
334 //
335 for (int i = 0; i != nImage; ++i) {
336 double const flux_i = fluxes[i];
337
338 for (int y = 0; y != height; ++y) {
339 typename ImageT::x_iterator iptr = imageList[i]->row_begin(y);
340 for (typename ImageT::Image::x_iterator mptr = mean.row_begin(y), end = mean.row_end(y);
341 mptr != end; ++mptr, ++iptr) {
342 if ((iptr.mask() & mask)) {
343 double const delta = ::fabs(flux_i * (*mptr) - iptr.image());
344 if (delta > maxChange) {
345 maxChange = delta;
346 }
347 iptr.image() = flux_i * (*mptr);
348 }
349 }
350 }
351 }
352 } else {
353 if (ncomp > static_cast<int>(eigenImages.size())) {
355 (boost::format("You only have %d eigen images (you asked for %d)") %
356 eigenImages.size() % ncomp)
357 .str());
358 }
359
360 for (int i = 0; i != nImage; ++i) {
362 fitEigenImagesToImage(eigenImages, ncomp, *imageList[i]);
363
364 for (int y = 0; y != height; ++y) {
365 typename ImageT::x_iterator iptr = imageList[i]->row_begin(y);
366 for (typename ImageT::Image::const_x_iterator fptr = fitted->row_begin(y),
367 end = fitted->row_end(y);
368 fptr != end; ++fptr, ++iptr) {
369 if (iptr.mask() & mask) {
370 double const delta = fabs(static_cast<double>(*fptr) - iptr.image());
371 if (delta > maxChange) {
372 maxChange = delta;
373 }
374
375 iptr.image() = *fptr;
376 }
377 }
378 }
379 }
380 }
381
382 return maxChange;
383}
384} // namespace
385template <typename ImageT>
386double ImagePca<ImageT>::updateBadPixels(unsigned long mask, int const ncomp) {
387 return do_updateBadPixels<ImageT>(typename ImageT::image_category(), _imageList, _fluxList, _eigenImages,
388 mask, ncomp);
389}
390
391namespace {
392template <typename T, typename U>
393struct IsSame {
394 IsSame(T const&, U const&) {}
395 bool operator()() { return false; }
396};
397
398template <typename T>
399struct IsSame<T, T> {
400 IsSame(T const& im1, T const& im2) : _same(im1.row_begin(0) == im2.row_begin(0)) {}
401 bool operator()() { return _same; }
402
403private:
404 bool _same;
405};
406
407// Test if two Images are identical; they need not be of the same type
408template <typename Image1T, typename Image2T>
409bool imagesAreIdentical(Image1T const& im1, Image2T const& im2) {
410 return IsSame<Image1T, Image2T>(im1, im2)();
411}
412} // namespace
413template <typename Image1T, typename Image2T>
414double innerProduct(Image1T const& lhs, Image2T const& rhs, int border) {
415 if (lhs.getWidth() <= 2 * border || lhs.getHeight() <= 2 * border) {
417 (boost::format("All image pixels are in the border of width %d: %dx%d") % border %
418 lhs.getWidth() % lhs.getHeight())
419 .str());
420 }
421
422 double sum = 0.0;
423 //
424 // Handle I.I specially for efficiency, and to avoid advancing the iterator twice
425 //
426 if (imagesAreIdentical(lhs, rhs)) {
427 for (int y = border; y != lhs.getHeight() - border; ++y) {
428 for (typename Image1T::const_x_iterator lptr = lhs.row_begin(y) + border,
429 lend = lhs.row_end(y) - border;
430 lptr != lend; ++lptr) {
431 typename Image1T::Pixel val = *lptr;
432 if (std::isfinite(val)) {
433 sum += val * val;
434 }
435 }
436 }
437 } else {
438 if (lhs.getDimensions() != rhs.getDimensions()) {
440 (boost::format("Dimension mismatch: %dx%d v. %dx%d") % lhs.getWidth() %
441 lhs.getHeight() % rhs.getWidth() % rhs.getHeight())
442 .str());
443 }
444
445 for (int y = border; y != lhs.getHeight() - border; ++y) {
446 typename Image2T::const_x_iterator rptr = rhs.row_begin(y) + border;
447 for (typename Image1T::const_x_iterator lptr = lhs.row_begin(y) + border,
448 lend = lhs.row_end(y) - border;
449 lptr != lend; ++lptr, ++rptr) {
450 double const tmp = (*lptr) * (*rptr);
451 if (std::isfinite(tmp)) {
452 sum += tmp;
453 }
454 }
455 }
456 }
457
458 return sum;
459}
460
461//
462// Explicit instantiations
463//
465#define INSTANTIATE(T) \
466 template class ImagePca<Image<T> >; \
467 template double innerProduct(Image<T> const&, Image<T> const&, int); \
468 template class ImagePca<MaskedImage<T> >;
469
470#define INSTANTIATE2(T, U) \
471 template double innerProduct(Image<T> const&, Image<U> const&, int); \
472 template double innerProduct(Image<U> const&, Image<T> const&, int);
473
476INSTANTIATE(int)
477INSTANTIATE(float)
478INSTANTIATE(double)
479
480INSTANTIATE2(float, double) // the two types must be different
482} // namespace image
483} // namespace afw
484} // namespace lsst
int end
double x
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
afw::table::PointKey< int > dimensions
Definition: GaussianPsf.cc:48
afw::table::Key< afw::table::Array< MaskPixelT > > mask
#define INSTANTIATE2(ImagePixelT1, ImagePixelT2)
Definition: MaskedImage.cc:690
uint64_t * ptr
Definition: RangeSet.cc:88
int y
Definition: SpanSet.cc:48
table::Key< int > b
table::Key< int > a
T begin(T... args)
A class to represent a 2-dimensional array of pixels.
Definition: Image.h:51
virtual double updateBadPixels(unsigned long mask, int const ncomp)
Update the bad pixels (i.e.
Definition: ImagePca.cc:386
void addImage(std::shared_ptr< ImageT > img, double flux=0.0)
Add an image to the set to be analyzed.
Definition: ImagePca.cc:64
ImageList getImageList() const
Return the list of images being analyzed.
Definition: ImagePca.cc:85
ImagePca(bool constantWeight=true)
ctor
Definition: ImagePca.cc:46
std::shared_ptr< ImageT > getMean() const
Return the mean of the images in ImagePca's list.
Definition: ImagePca.cc:90
std::vector< std::shared_ptr< ImageT > > ImageList
Definition: ImagePca.h:47
ImagePca & operator=(ImagePca const &)
virtual void analyze()
Definition: ImagePca.cc:123
Reports attempts to exceed implementation-defined length limits for some classes.
Definition: Runtime.h:76
Reports attempts to access elements outside a valid range of indices.
Definition: Runtime.h:89
T emplace_back(T... args)
T end(T... args)
T fabs(T... args)
T isfinite(T... args)
#define INSTANTIATE(TYPE)
Definition: ImagePca.cc:127
def dot(symb, c, r, frame=None, size=2, ctype=None, origin=afwImage.PARENT, *args, **kwargs)
Definition: ds9.py:100
Backwards-compatibility support for depersisting the old Calib (FluxMag0/FluxMag0Err) objects.
double innerProduct(Image1T const &lhs, Image2T const &rhs, int const border=0)
Calculate the inner product of two images.
Definition: ImagePca.cc:414
float Pixel
Typedefs to be used for pixel values.
Definition: common.h:37
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174
A base class for image defects.
STL namespace.
T reserve(T... args)
T size(T... args)
T sort(T... args)
afw::table::Key< double > weight
ImageT val
Definition: CR.cc:146