LSSTApplications  10.0+286,10.0+36,10.0+46,10.0-2-g4f67435,10.1+152,10.1+37,11.0,11.0+1,11.0-1-g47edd16,11.0-1-g60db491,11.0-1-g7418c06,11.0-2-g04d2804,11.0-2-g68503cd,11.0-2-g818369d,11.0-2-gb8b8ce7
LSSTDataManagementBasePackage
Public Types | Public Member Functions | Private Attributes | List of all members
lsst::meas::algorithms::Shapelet Class Reference

#include <Shapelet.h>

Public Types

typedef float PixelT
 This class includes the basic functionality for shapelets. More...
 
typedef boost::shared_ptr
< Shapelet
Ptr
 
typedef boost::shared_ptr
< const Shapelet
ConstPtr
 
typedef Eigen::VectorXd ShapeletVector
 
typedef Eigen::MatrixXd ShapeletCovariance
 
typedef
lsst::afw::table::SourceRecord 
Source
 
typedef
lsst::afw::image::Exposure
< PixelT
Exposure
 
typedef lsst::afw::geom::PointD PointD
 
typedef lsst::afw::image::MaskPixel MaskPixel
 

Public Member Functions

 Shapelet (int order, double sigma)
 Basic constructor requires order and sigma. More...
 
 Shapelet (int order, double sigma, const ShapeletVector &vector)
 A constructor from a vector of values. More...
 
 Shapelet (int order, double sigma, const ShapeletVector &vector, const ShapeletCovariance &cov)
 A constructor from a vector of values, with covariance. More...
 
 ~Shapelet ()
 Destructor needs to delete pImpl. More...
 
 Shapelet (const Shapelet &rhs)
 Copy constructor does a deep copy. More...
 
Shapeletoperator= (const Shapelet &rhs)
 op= does a deep copy More...
 
int getOrder () const
 get the order of the shapelet More...
 
double getSigma () const
 get the scale size of the shapelet More...
 
int size () const
 the size of the shapelet vector More...
 
const ShapeletVectorgetValues () const
 get the values as a vector. More...
 
bool hasCovariance () const
 does the shapelet have a covariance matrix stored? More...
 
boost::shared_ptr< const
ShapeletCovariance
getCovariance () const
 get the covariance matrix More...
 
void setSigma (double sigma)
 set a new value of sigma More...
 
std::complex< double > getPQ (int p, int q)
 Get a complex b_pq value. More...
 
double evaluateAt (const PointD &pos)
 Evaluate f(x,y) More...
 
double evaluateAt (double x, double y)
 
bool measureFromImage (const Source &source, const PointD &pos, bool isCentroidFixed, bool isSigmaFixed, double aperture, const Exposure &exposure, const MaskPixel okmask=0)
 measure shapelet decomposition of an image More...
 
 Shapelet (const shapelet::BVec &bvec)
 A constructor that is only used by various Shapelet implementations. More...
 
const shapelet::BVecviewAsBVec () const
 View the shapelet as a BVec. More...
 
shapelet::BVecviewAsBVec ()
 

Private Attributes

ShapeletImplpImpl
 

Detailed Description

Definition at line 50 of file Shapelet.h.

Member Typedef Documentation

typedef boost::shared_ptr<const Shapelet> lsst::meas::algorithms::Shapelet::ConstPtr

Definition at line 95 of file Shapelet.h.

Definition at line 101 of file Shapelet.h.

Definition at line 103 of file Shapelet.h.

This class includes the basic functionality for shapelets.

There are a few definitions of shapelets out there. This class implements the shapelets defined in Bernstein & Jarvis (2002) (where they are called the Gauss-Laguerre basis functions). These are similar to the "polar shapelets" of Massey & Refregier.

Shapelets are a complete basis set that can be used to describe any image. However, they are most useful for describing things that are intrinsically similar to a 2-d Gaussian. Since the PSF (for ground-based telescopes) is dominated by a Gaussian component, this makes them a good choice for describing the PSF.

For specificity, the basis functions we use are:

For p >= q: psi_pq(x,y,sigma) = (pi p! q!)^-1/2 z^m exp(-r^2/2) K_pq(r^2) where z = (x + Iy)/sigma m = p-q r = |z|^2 K_pq(r^2) = (-)^q q! L_q^(m) (r^2) L_q^(m)(x) are Laguerre polynomials

For any particular value of sigma, these define a complete set of basis functions.

The image is thus decomposed as:

f(x,y) = Sum_pq b_pq psi_pq(x,y)

As with any basis set, the decomposition is only exact if an infinite number of basis functions are used. We truncate the expansion at p+q <= order for some (specified) order.

Definition at line 93 of file Shapelet.h.

Definition at line 102 of file Shapelet.h.

typedef boost::shared_ptr<Shapelet> lsst::meas::algorithms::Shapelet::Ptr

Definition at line 94 of file Shapelet.h.

Definition at line 98 of file Shapelet.h.

Definition at line 97 of file Shapelet.h.

Definition at line 100 of file Shapelet.h.

Constructor & Destructor Documentation

lsst::meas::algorithms::Shapelet::Shapelet ( int  order,
double  sigma 
)

Basic constructor requires order and sigma.

order defines how many shapelet coefficients to measure. e.g. order = 4 includes: b00, b10, b20, b11, b30, b21, b40, b31, b22 (bii are real, the others are complex.)

sigma is the scale size of Gaussian factor in the shapelet decomposition, measured in arcsec.

Parameters
orderthe maximum value of p+q to measure.
sigmathe scale size of the shapelet decomposition (arcsec).

Definition at line 145 of file Shapelet.cc.

145  :
146  pImpl(new ShapeletImpl(order,sigma))
147  {}
afw::table::Key< double > sigma
Definition: GaussianPsf.cc:43
lsst::meas::algorithms::Shapelet::Shapelet ( int  order,
double  sigma,
const ShapeletVector vector 
)

A constructor from a vector of values.

The input vector should have (order+1)*(order+2)/2 elements.

Caveat: The input vector is real, not complex. See the comment for getValues for the expected order of values.

Parameters
orderthe maximum value of p+q to measure.
sigmathe scale size of the shapelet decomposition (arcsec).
vectorthe shapelet vector

Definition at line 149 of file Shapelet.cc.

149  :
150  pImpl(new ShapeletImpl(order,sigma,vector))
151  {}
afw::table::Key< double > sigma
Definition: GaussianPsf.cc:43
lsst::meas::algorithms::Shapelet::Shapelet ( int  order,
double  sigma,
const ShapeletVector vector,
const ShapeletCovariance cov 
)

A constructor from a vector of values, with covariance.

The input vector should have (order+1)*(order+2)/2 elements. The covariance matrix is of the input vector of values. So it should be square, symmetric, and have the same size in each dimension as the input vector.

Caveat: The input vector is real, not complex. See the comment for getValues for the expected order of values.

Parameters
orderthe maximum value of p+q to measure.
sigmathe scale size of the shapelet decomposition (arcsec).
vectorthe shapelet vector
covthe covariance matrix

Definition at line 153 of file Shapelet.cc.

153  :
154  pImpl(new ShapeletImpl(order,sigma,vector,cov))
155  {}
afw::table::Key< double > sigma
Definition: GaussianPsf.cc:43
lsst::meas::algorithms::Shapelet::~Shapelet ( )

Destructor needs to delete pImpl.

Definition at line 161 of file Shapelet.cc.

162  { delete pImpl; pImpl = 0; }
lsst::meas::algorithms::Shapelet::Shapelet ( const Shapelet rhs)

Copy constructor does a deep copy.

Definition at line 157 of file Shapelet.cc.

157  :
158  pImpl(new ShapeletImpl(*rhs.pImpl))
159  {}
lsst::meas::algorithms::Shapelet::Shapelet ( const shapelet::BVec bvec)

A constructor that is only used by various Shapelet implementations.

This is used by ShapeletInterpolation for example. It is more convenient to leave this public, but you can't actually use it, since BVec isn't defined publicly.

Definition at line 408 of file Shapelet.cc.

408 : pImpl(new ShapeletImpl(bvec)) {}

Member Function Documentation

double lsst::meas::algorithms::Shapelet::evaluateAt ( const PointD pos)

Evaluate f(x,y)

Evaluate the intensity as a function of chip position (x,y)

This is the approximation to the true f(x,y) defined by the decomposition Sum_pq b_pq psi_pq(x,y,sigma).

Definition at line 194 of file Shapelet.cc.

195  { return evaluateAt(pos.getX(),pos.getY()); }
double evaluateAt(const PointD &pos)
Evaluate f(x,y)
Definition: Shapelet.cc:194
double lsst::meas::algorithms::Shapelet::evaluateAt ( double  x,
double  y 
)

Definition at line 197 of file Shapelet.cc.

198  {
199  // TODO: MJ - I couldn't reproduce the crash with g++ 4.3.3.
200  // Need to try this on lsst computers I guess.
201 #if 0 // crashes g++ 4.3.5
202  return pImpl->evaluateAt(x,y);
203 #else
204  // MJ - This is definitely wrong! It starts an infinite loop!
206 #endif
207  }
int y
double evaluateAt(const PointD &pos)
Evaluate f(x,y)
Definition: Shapelet.cc:194
double evaluateAt(double x, double y)
Definition: Shapelet.cc:114
double x
boost::shared_ptr< const Shapelet::ShapeletCovariance > lsst::meas::algorithms::Shapelet::getCovariance ( ) const

get the covariance matrix

Definition at line 185 of file Shapelet.cc.

186  { return pImpl->getCovariance(); }
boost::shared_ptr< ShapeletCovariance > & getCovariance()
Definition: Shapelet.cc:94
int lsst::meas::algorithms::Shapelet::getOrder ( ) const

get the order of the shapelet

Definition at line 170 of file Shapelet.cc.

171  { return pImpl->getOrder(); }
std::complex< double > lsst::meas::algorithms::Shapelet::getPQ ( int  p,
int  q 
)

Get a complex b_pq value.

Get the coefficient b_pq for the shapelet decomposition.

Definition at line 191 of file Shapelet.cc.

192  { return pImpl->getPQ(p,q); }
std::complex< double > getPQ(int p, int q)
Definition: Shapelet.cc:97
double lsst::meas::algorithms::Shapelet::getSigma ( ) const

get the scale size of the shapelet

Definition at line 173 of file Shapelet.cc.

174  { return pImpl->getSigma(); }
const Shapelet::ShapeletVector & lsst::meas::algorithms::Shapelet::getValues ( ) const

get the values as a vector.

The order of values are: b00 real(b10) imag(b10) real(b20) imag(b20) b11 real(b30) imag(b30) real(b21) imag(b21) real(b40) imag(b40) real(b31) imag(b31) b22 etc.

where bpq corresponds to the coefficients $b_{pq}$ as defined in Bernstein and Jarvis (2002).

Definition at line 179 of file Shapelet.cc.

180  { return pImpl->getValues(); }
const DVector & getValues() const
Definition: BVec.h:66
bool lsst::meas::algorithms::Shapelet::hasCovariance ( ) const

does the shapelet have a covariance matrix stored?

Definition at line 182 of file Shapelet.cc.

183  { return pImpl->getCovariance().get(); }
boost::shared_ptr< ShapeletCovariance > & getCovariance()
Definition: Shapelet.cc:94
bool lsst::meas::algorithms::Shapelet::measureFromImage ( const Source source,
const PointD pos,
bool  isCentroidFixed,
bool  isSigmaFixed,
double  aperture,
const Exposure exposure,
const MaskPixel  okmask = 0 
)

measure shapelet decomposition of an image

The initial estimate of the centroid is taken from the source. The initial estimate of sigma is the value given by this->getSigma().

The normal operation would be to allow both values to shift to find the best values for the image. However, there are reasons to keep each fixed.

First, if we have very accurate positions given in the source object, and we want the PSF to include the centroid shift that is present in the seeing, then this would be a good reason to keep it fixed.

Second, interpolation of the PSF is more straightforward if all the stars are measured with the saem sigma. So the normal practice is to let sigma float for each star, then find the mean of all the stars in an image, and finally remeasure each star with this constant value of sigma.

So, there are two boolean parameters that govern what should be kept fixed during the calculation.

isCentroidFixed says whether to keep the centoid fixed. If it is false, then the centoid is allowed to vary until b10 = 0 (our definition of being centroided).

isSigmaFixed says whether to keep the sigma fixed. If it is false, then sigma is allowed to vary until b11 = 0 (which corresponds to when the measured shapelet components have the maximum signal-to-noise).

The measurement is done using a circular aperture, whose radius is given by the paremeter aperture. This aperture (given in arcsec) defines a circle in sky coordinates, not chip coordinates.

The return value is true if the measurement is successful, and falso if not.

Parameters
sourceThe source of which to measure the decomposition.
posThe position of the object
isCentroidFixedIs sigma fixed? or should it be allowed to vary?
isSigmaFixedIs sigma fixed? or should it be allowed to vary?
apertureThe aperture size in arcsec for the measurement.
exposureThe exposure on which to measure the decompoisition.
okmaskThe mask values that are ok to use

Definition at line 358 of file Shapelet.cc.

363  {
364  using shapelet::Ellipse;
365  using shapelet::PixelList;
366 
367  std::vector<PixelList> pix(1);
368  // Fill PixelList with pixel data around position pos:
369  getPixList(pix[0], source, pos, aperture, exposure, okmask);
370 
371  double sigma = pImpl->getSigma();
372  Ellipse ell;
373  ell.fixGam();
374  if (isCentroidFixed) ell.fixCen();
375  else {
376  // Initial crude estimates to get close to the right value
377  // in case we have a poor starting point.
378  // TODO: These might not be necessary for LSST.
379  // Should compare speed with and without this step.
380  ell.peakCentroid(pix[0],aperture/3.);
381  ell.crudeMeasure(pix[0],sigma);
382  }
383  if (isSigmaFixed) ell.fixMu();
384  long flag = 0;
385  int order = getOrder();
386  if (!isCentroidFixed || !isSigmaFixed) {
387  if (!ell.measure(pix,order,order+4,order,sigma,flag,1.e-4)) {
388  return false;
389  }
390  if (flag) return false;
391  if (!isSigmaFixed) {
392  double mu = ell.getMu();
393  sigma *= exp(mu);
394  dbg<<"sigma = "<<sigma<<std::endl;
395  Assert(sigma > 0);
396  pImpl->setSigma(sigma);
397  }
398  }
399 
400  ShapeletCovariance* cov = pImpl->getCovariance().get();
401  if (ell.measureShapelet(pix,*pImpl,order,order+4,order,cov)) {
402  return true;
403  } else {
404  return false;
405  }
406  }
int getOrder() const
get the order of the shapelet
Definition: Shapelet.cc:170
Eigen::MatrixXd ShapeletCovariance
Definition: Shapelet.h:98
afw::table::Key< double > sigma
Definition: GaussianPsf.cc:43
void setSigma(double sigma)
Definition: BVec.h:69
#define dbg
Definition: dbg.h:69
boost::shared_ptr< ShapeletCovariance > & getCovariance()
Definition: Shapelet.cc:94
#define Assert(x)
Definition: dbg.h:73
Shapelet & lsst::meas::algorithms::Shapelet::operator= ( const Shapelet rhs)

op= does a deep copy

Definition at line 164 of file Shapelet.cc.

165  {
166  *pImpl = *rhs.pImpl;
167  return *this;
168  }
void lsst::meas::algorithms::Shapelet::setSigma ( double  sigma)

set a new value of sigma

Definition at line 188 of file Shapelet.cc.

189  { return pImpl->setSigma(sigma); }
afw::table::Key< double > sigma
Definition: GaussianPsf.cc:43
void setSigma(double sigma)
Definition: BVec.h:69
int lsst::meas::algorithms::Shapelet::size ( ) const

the size of the shapelet vector

Definition at line 176 of file Shapelet.cc.

177  { return pImpl->size(); }
const shapelet::BVec & lsst::meas::algorithms::Shapelet::viewAsBVec ( ) const

View the shapelet as a BVec.

This is also only used by various Shapelet implementations, but it is convenient to leave it public rather than define a bunch of friend classes.

Definition at line 410 of file Shapelet.cc.

410 { return *pImpl; }
shapelet::BVec & lsst::meas::algorithms::Shapelet::viewAsBVec ( )

Definition at line 411 of file Shapelet.cc.

411 { return *pImpl; }

Member Data Documentation

ShapeletImpl* lsst::meas::algorithms::Shapelet::pImpl
private

Definition at line 306 of file Shapelet.h.


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