LSSTApplications  10.0-2-g4f67435,11.0.rc2+1,11.0.rc2+12,11.0.rc2+3,11.0.rc2+4,11.0.rc2+5,11.0.rc2+6,11.0.rc2+7,11.0.rc2+8
LSSTDataManagementBasePackage
Public Types | Public Member Functions | Private Member Functions | Private Attributes | 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

typedef boost::shared_ptr< ImageT > Ptr
 
typedef boost::shared_ptr
< const ImageT > 
ConstPtr
 
typedef std::vector< typename
ImageT::Ptr > 
ImageList
 

Public Member Functions

 ImagePca (bool constantWeight=true)
 ctor More...
 
virtual ~ImagePca ()
 
void addImage (typename ImageT::Ptr img, double flux=0.0)
 
ImageList getImageList () const
 Return the list of images being analyzed. More...
 
geom::Extent2I const getDimensions () const
 Return the dimension of the images being analyzed. More...
 
ImageT::Ptr getMean () const
 
virtual void analyze ()
 
virtual double updateBadPixels (unsigned long mask, int const ncomp)
 
std::vector< double > const & getEigenValues () const
 Return Eigen values. More...
 
ImageList const & getEigenImages () const
 Return Eigen images. More...
 

Private Member Functions

double getFlux (int i) const
 

Private Attributes

ImageList _imageList
 
std::vector< double > _fluxList
 
geom::Extent2I _dimensions
 
bool _constantWeight
 
std::vector< double > _eigenValues
 
ImageList _eigenImages
 

Detailed Description

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

Definition at line 47 of file ImagePca.h.

Member Typedef Documentation

template<typename ImageT>
typedef boost::shared_ptr<const ImageT> lsst.afw.image::ImagePca< ImageT >::ConstPtr

Definition at line 50 of file ImagePca.h.

template<typename ImageT>
typedef std::vector<typename ImageT::Ptr> lsst.afw.image::ImagePca< ImageT >::ImageList

Definition at line 52 of file ImagePca.h.

template<typename ImageT>
typedef boost::shared_ptr<ImageT> lsst.afw.image::ImagePca< ImageT >::Ptr

Definition at line 49 of file ImagePca.h.

Constructor & Destructor Documentation

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

ctor

Parameters
constantWeightShould all stars be weighted equally?

Definition at line 49 of file ImagePca.cc.

50  :
51  _imageList(),
52  _fluxList(),
53  _dimensions(0,0),
54  _constantWeight(constantWeight),
55  _eigenValues(std::vector<double>()),
57 }
ImageList _eigenImages
Definition: ImagePca.h:82
std::vector< double > _fluxList
Definition: ImagePca.h:76
std::vector< double > _eigenValues
Definition: ImagePca.h:81
geom::Extent2I _dimensions
Definition: ImagePca.h:77
std::vector< typename ImageT::Ptr > ImageList
Definition: ImagePca.h:52
template<typename ImageT>
virtual lsst.afw.image::ImagePca< ImageT >::~ImagePca ( )
inlinevirtual

Definition at line 55 of file ImagePca.h.

55 {}

Member Function Documentation

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

Add an image to the set to be analyzed

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

Definition at line 65 of file ImagePca.cc.

67  {
68  if (_imageList.empty()) {
69  _dimensions = img->getDimensions();
70  } else {
71  if (getDimensions() != img->getDimensions()) {
72  throw LSST_EXCEPT(
73  lsst::pex::exceptions::LengthError,
74  (boost::format("Dimension mismatch: saw %dx%d; expected %dx%d") %
75  img->getWidth() % img->getHeight() %
76  _dimensions.getX() % _dimensions.getY()
77  ).str()
78  );
79  }
80  }
81 
82  if (flux == 0.0) {
83  throw LSST_EXCEPT(lsst::pex::exceptions::OutOfRangeError, "Flux may not be zero");
84  }
85 
86  _imageList.push_back(img);
87  _fluxList.push_back(flux);
88 }
std::vector< double > _fluxList
Definition: ImagePca.h:76
geom::Extent2I _dimensions
Definition: ImagePca.h:77
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46
geom::Extent2I const getDimensions() const
Return the dimension of the images being analyzed.
Definition: ImagePca.h:61
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 135 of file ImagePca.cc.

136 {
137  int const nImage = _imageList.size();
138 
139  if (nImage == 0) {
140  throw LSST_EXCEPT(lsst::pex::exceptions::LengthError, "No images provided for PCA analysis");
141  }
142  /*
143  * Eigen doesn't like 1x1 matrices, but we don't really need it to handle a single matrix...
144  */
145  if (nImage == 1) {
146  _eigenImages.clear();
147  _eigenImages.push_back(typename ImageT::Ptr(new ImageT(*_imageList[0], true)));
148 
149  _eigenValues.clear();
150  _eigenValues.push_back(1.0);
151 
152  return;
153  }
154  /*
155  * Find the eigenvectors/values of the scalar product matrix, R' (Eq. 7.4)
156  */
157  Eigen::MatrixXd R(nImage, nImage); // residuals' inner products
158 
159  double flux_bar = 0; // mean of flux for all regions
160  for (int i = 0; i != nImage; ++i) {
161  typename GetImage<ImageT>::type const& im_i = *GetImage<ImageT>::getImage(_imageList[i]);
162  double const flux_i = getFlux(i);
163  flux_bar += flux_i;
164 
165  for (int j = i; j != nImage; ++j) {
166  typename GetImage<ImageT>::type const& im_j = *GetImage<ImageT>::getImage(_imageList[j]);
167  double const flux_j = getFlux(j);
168 
169  double dot = innerProduct(im_i, im_j);
170  if (_constantWeight) {
171  dot /= flux_i*flux_j;
172  }
173  R(i, j) = R(j, i) = dot/nImage;
174  }
175  }
176  flux_bar /= nImage;
177  Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eVecValues(R);
178  Eigen::MatrixXd const& Q = eVecValues.eigenvectors();
179  Eigen::VectorXd const& lambda = eVecValues.eigenvalues();
180  //
181  // We need to sort the eigenValues, and remember the permutation we applied to the eigenImages
182  // We'll use the vector lambdaAndIndex to achieve this
183  //
184  std::vector<std::pair<double, int> > lambdaAndIndex; // pairs (eValue, index)
185  lambdaAndIndex.reserve(nImage);
186 
187  for (int i = 0; i != nImage; ++i) {
188  lambdaAndIndex.push_back(std::make_pair(lambda(i), i));
189  }
190  std::sort(lambdaAndIndex.begin(), lambdaAndIndex.end(), SortEvalueDecreasing<double>());
191  //
192  // Save the (sorted) eigen values
193  //
194  _eigenValues.clear();
195  _eigenValues.reserve(nImage);
196  for (int i = 0; i != nImage; ++i) {
197  _eigenValues.push_back(lambdaAndIndex[i].first);
198  }
199  //
200  // Contruct the first ncomp eigenimages in basis
201  //
202  int ncomp = 100; // number of components to keep
203 
204  _eigenImages.clear();
205  _eigenImages.reserve(ncomp < nImage ? ncomp : nImage);
206 
207  for(int i = 0; i < ncomp; ++i) {
208  if(i >= nImage) {
209  continue;
210  }
211 
212  int const ii = lambdaAndIndex[i].second; // the index after sorting (backwards) by eigenvalue
213 
214  typename ImageT::Ptr eImage(new ImageT(_dimensions));
215  *eImage = static_cast<typename ImageT::Pixel>(0);
216 
217  for (int j = 0; j != nImage; ++j) {
218  int const jj = lambdaAndIndex[j].second; // the index after sorting (backwards) by eigenvalue
219  double const weight = Q(jj, ii)*(_constantWeight ? flux_bar/getFlux(jj) : 1);
220  eImage->scaledPlus(weight, *_imageList[jj]);
221  }
222  _eigenImages.push_back(eImage);
223  }
224 }
tbl::Key< double > weight
double getFlux(int i) const
Definition: ImagePca.h:73
double innerProduct(Image1T const &lhs, Image2T const &rhs, int const border=0)
Definition: ImagePca.cc:454
ImageList _eigenImages
Definition: ImagePca.h:82
std::vector< double > _eigenValues
Definition: ImagePca.h:81
geom::Extent2I _dimensions
Definition: ImagePca.h:77
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46
template<typename ImageT>
geom::Extent2I const lsst.afw.image::ImagePca< ImageT >::getDimensions ( ) const
inline

Return the dimension of the images being analyzed.

Definition at line 61 of file ImagePca.h.

61 { return _dimensions; }
geom::Extent2I _dimensions
Definition: ImagePca.h:77
template<typename ImageT>
ImageList const& lsst.afw.image::ImagePca< ImageT >::getEigenImages ( ) const
inline

Return Eigen images.

Definition at line 70 of file ImagePca.h.

70 { return _eigenImages; }
ImageList _eigenImages
Definition: ImagePca.h:82
template<typename ImageT>
std::vector<double> const& lsst.afw.image::ImagePca< ImageT >::getEigenValues ( ) const
inline

Return Eigen values.

Definition at line 68 of file ImagePca.h.

68 { return _eigenValues; }
std::vector< double > _eigenValues
Definition: ImagePca.h:81
template<typename ImageT>
double lsst.afw.image::ImagePca< ImageT >::getFlux ( int  i) const
inlineprivate

Definition at line 73 of file ImagePca.h.

73 { return _fluxList[i]; }
std::vector< double > _fluxList
Definition: ImagePca.h:76
template<typename ImageT >
ImagePca< ImageT >::ImageList lsst.afw.image::ImagePca< ImageT >::getImageList ( ) const

Return the list of images being analyzed.

Definition at line 92 of file ImagePca.cc.

92  {
93  return _imageList;
94 }
template<typename ImageT >
ImageT::Ptr lsst.afw.image::ImagePca< ImageT >::getMean ( ) const

Return the mean of the images in ImagePca's list

Definition at line 101 of file ImagePca.cc.

101  {
102  if (_imageList.empty()) {
103  throw LSST_EXCEPT(lsst::pex::exceptions::LengthError, "You haven't provided any images");
104  }
105 
106  typename ImageT::Ptr mean(new ImageT(getDimensions()));
107  *mean = static_cast<typename ImageT::Pixel>(0);
108 
109  for (typename ImageList::const_iterator ptr = _imageList.begin(), end = _imageList.end();
110  ptr != end; ++ptr) {
111  *mean += **ptr;
112  }
113  *mean /= _imageList.size();
114 
115  return mean;
116 }
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46
geom::Extent2I const getDimensions() const
Return the dimension of the images being analyzed.
Definition: ImagePca.h:61
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

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

Parameters
maskMask defining bad pixels
ncompNumber of components to use in estimate

Definition at line 417 of file ImagePca.cc.

421 {
422  return do_updateBadPixels<ImageT>(typename ImageT::image_category(),
423  _imageList, _fluxList, _eigenImages, mask, ncomp);
424 }
ImageList _eigenImages
Definition: ImagePca.h:82
std::vector< double > _fluxList
Definition: ImagePca.h:76

Member Data Documentation

template<typename ImageT>
bool lsst.afw.image::ImagePca< ImageT >::_constantWeight
private

Definition at line 79 of file ImagePca.h.

template<typename ImageT>
geom::Extent2I lsst.afw.image::ImagePca< ImageT >::_dimensions
private

Definition at line 77 of file ImagePca.h.

template<typename ImageT>
ImageList lsst.afw.image::ImagePca< ImageT >::_eigenImages
private

Definition at line 82 of file ImagePca.h.

template<typename ImageT>
std::vector<double> lsst.afw.image::ImagePca< ImageT >::_eigenValues
private

Definition at line 81 of file ImagePca.h.

template<typename ImageT>
std::vector<double> lsst.afw.image::ImagePca< ImageT >::_fluxList
private

Definition at line 76 of file ImagePca.h.

template<typename ImageT>
ImageList lsst.afw.image::ImagePca< ImageT >::_imageList
private

Definition at line 75 of file ImagePca.h.


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