LSST Applications  21.0.0-172-gfb10e10a+18fedfabac,22.0.0+297cba6710,22.0.0+80564b0ff1,22.0.0+8d77f4f51a,22.0.0+a28f4c53b1,22.0.0+dcf3732eb2,22.0.1-1-g7d6de66+2a20fdde0d,22.0.1-1-g8e32f31+297cba6710,22.0.1-1-geca5380+7fa3b7d9b6,22.0.1-12-g44dc1dc+2a20fdde0d,22.0.1-15-g6a90155+515f58c32b,22.0.1-16-g9282f48+790f5f2caa,22.0.1-2-g92698f7+dcf3732eb2,22.0.1-2-ga9b0f51+7fa3b7d9b6,22.0.1-2-gd1925c9+bf4f0e694f,22.0.1-24-g1ad7a390+a9625a72a8,22.0.1-25-g5bf6245+3ad8ecd50b,22.0.1-25-gb120d7b+8b5510f75f,22.0.1-27-g97737f7+2a20fdde0d,22.0.1-32-gf62ce7b1+aa4237961e,22.0.1-4-g0b3f228+2a20fdde0d,22.0.1-4-g243d05b+871c1b8305,22.0.1-4-g3a563be+32dcf1063f,22.0.1-4-g44f2e3d+9e4ab0f4fa,22.0.1-42-gca6935d93+ba5e5ca3eb,22.0.1-5-g15c806e+85460ae5f3,22.0.1-5-g58711c4+611d128589,22.0.1-5-g75bb458+99c117b92f,22.0.1-6-g1c63a23+7fa3b7d9b6,22.0.1-6-g50866e6+84ff5a128b,22.0.1-6-g8d3140d+720564cf76,22.0.1-6-gd805d02+cc5644f571,22.0.1-8-ge5750ce+85460ae5f3,master-g6e05de7fdc+babf819c66,master-g99da0e417a+8d77f4f51a,w.2021.48
LSST Data Management Base Package
Public Types | Public Member Functions | List of all members
lsst::afw::image::ImagePca< ImageT > Class Template Reference

#include <ImagePca.h>

Inheritance diagram for lsst::afw::image::ImagePca< ImageT >:
lsst::ip::diffim::detail::KernelPca< ImageT > lsst::meas::algorithms::PsfImagePca< ImageT >

Public Types

using ImageList = std::vector< std::shared_ptr< ImageT > >
 

Public Member Functions

 ImagePca (bool constantWeight=true)
 ctor More...
 
virtual ~ImagePca ()=default
 
 ImagePca (ImagePca const &)
 
 ImagePca (ImagePca &&)
 
ImagePcaoperator= (ImagePca const &)
 
ImagePcaoperator= (ImagePca &&)
 
void addImage (std::shared_ptr< ImageT > img, double flux=0.0)
 Add an image to the set to be analyzed. More...
 
ImageList getImageList () const
 Return the list of images being analyzed. More...
 
lsst::geom::Extent2I const getDimensions () const
 Return the dimension of the images being analyzed. More...
 
std::shared_ptr< ImageT > getMean () const
 Return the mean of the images in ImagePca's list. More...
 
virtual void analyze ()
 
virtual double updateBadPixels (unsigned long mask, int const ncomp)
 Update the bad pixels (i.e. More...
 
std::vector< double > const & getEigenValues () const
 Return Eigen values. More...
 
ImageList const & getEigenImages () const
 Return Eigen images. More...
 

Detailed Description

template<typename ImageT>
class lsst::afw::image::ImagePca< ImageT >

Definition at line 45 of file ImagePca.h.

Member Typedef Documentation

◆ ImageList

template<typename ImageT >
using lsst::afw::image::ImagePca< ImageT >::ImageList = std::vector<std::shared_ptr<ImageT> >

Definition at line 47 of file ImagePca.h.

Constructor & Destructor Documentation

◆ ImagePca() [1/3]

template<typename ImageT >
lsst::afw::image::ImagePca< ImageT >::ImagePca ( bool  constantWeight = true)
explicit

ctor

Parameters
constantWeightShould all stars be weighted equally?

Definition at line 46 of file ImagePca.cc.

47  : _imageList(),
48  _fluxList(),
49  _dimensions(0, 0),
50  _constantWeight(constantWeight),
51  _eigenValues(std::vector<double>()),
52  _eigenImages(ImageList()) {}
std::vector< std::shared_ptr< ImageT > > ImageList
Definition: ImagePca.h:47

◆ ~ImagePca()

template<typename ImageT >
virtual lsst::afw::image::ImagePca< ImageT >::~ImagePca ( )
virtualdefault

◆ ImagePca() [2/3]

template<typename ImageT >
lsst::afw::image::ImagePca< ImageT >::ImagePca ( ImagePca< ImageT > const &  )
default

◆ ImagePca() [3/3]

template<typename ImageT >
lsst::afw::image::ImagePca< ImageT >::ImagePca ( ImagePca< ImageT > &&  )
default

Member Function Documentation

◆ addImage()

template<typename ImageT >
void lsst::afw::image::ImagePca< ImageT >::addImage ( std::shared_ptr< ImageT >  img,
double  flux = 0.0 
)

Add an image to the set to be analyzed.

Parameters
imgImage to add to set
fluxImage's flux
Exceptions
lsst::pex::exceptions::LengthErrorif all the images aren't the same size

Definition at line 64 of file ImagePca.cc.

64  {
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 }
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
lsst::geom::Extent2I const getDimensions() const
Return the dimension of the images being analyzed.
Definition: ImagePca.h:75
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 empty(T... args)
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174
T push_back(T... args)

◆ analyze()

template<typename ImageT >
void lsst::afw::image::ImagePca< ImageT >::analyze
virtual

Reimplemented in lsst::meas::algorithms::PsfImagePca< ImageT >, and lsst::ip::diffim::detail::KernelPca< ImageT >.

Definition at line 123 of file ImagePca.cc.

123  {
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 }
T begin(T... args)
T clear(T... args)
T emplace_back(T... args)
T end(T... args)
def dot(symb, c, r, frame=None, size=2, ctype=None, origin=afwImage.PARENT, *args, **kwargs)
Definition: ds9.py:100
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
T reserve(T... args)
T size(T... args)
T sort(T... args)
afw::table::Key< double > weight

◆ getDimensions()

template<typename ImageT >
lsst::geom::Extent2I const lsst::afw::image::ImagePca< ImageT >::getDimensions ( ) const
inline

Return the dimension of the images being analyzed.

Definition at line 75 of file ImagePca.h.

75 { return _dimensions; }

◆ getEigenImages()

template<typename ImageT >
ImageList const& lsst::afw::image::ImagePca< ImageT >::getEigenImages ( ) const
inline

Return Eigen images.

Definition at line 100 of file ImagePca.h.

100 { return _eigenImages; }

◆ getEigenValues()

template<typename ImageT >
std::vector<double> const& lsst::afw::image::ImagePca< ImageT >::getEigenValues ( ) const
inline

Return Eigen values.

Definition at line 98 of file ImagePca.h.

98 { return _eigenValues; }

◆ getImageList()

template<typename ImageT >
ImagePca< ImageT >::ImageList lsst::afw::image::ImagePca< ImageT >::getImageList

Return the list of images being analyzed.

Definition at line 85 of file ImagePca.cc.

85  {
86  return _imageList;
87 }

◆ getMean()

template<typename ImageT >
std::shared_ptr< ImageT > lsst::afw::image::ImagePca< ImageT >::getMean

Return the mean of the images in ImagePca's list.

Definition at line 90 of file ImagePca.cc.

90  {
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 }
int end
uint64_t * ptr
Definition: RangeSet.cc:88

◆ operator=() [1/2]

template<typename ImageT >
ImagePca< ImageT > & lsst::afw::image::ImagePca< ImageT >::operator= ( ImagePca< ImageT > &&  )
default

◆ operator=() [2/2]

template<typename ImageT >
ImagePca< ImageT > & lsst::afw::image::ImagePca< ImageT >::operator= ( ImagePca< ImageT > const &  )
default

◆ updateBadPixels()

template<typename ImageT >
double lsst::afw::image::ImagePca< ImageT >::updateBadPixels ( unsigned long  mask,
int const  ncomp 
)
virtual

Update the bad pixels (i.e.

those for which (value & mask) != 0) based on the current PCA decomposition; if none is available, use the mean of the good pixels

Parameters
maskMask defining bad pixels
ncompNumber of components to use in estimate
Returns
the maximum change made to any pixel

N.b. the work is actually done in do_updateBadPixels as the code only makes sense and compiles when we are doing a PCA on a set of MaskedImages

Definition at line 386 of file ImagePca.cc.

386  {
387  return do_updateBadPixels<ImageT>(typename ImageT::image_category(), _imageList, _fluxList, _eigenImages,
388  mask, ncomp);
389 }
afw::table::Key< afw::table::Array< MaskPixelT > > mask

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