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
|
#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... | |
Shapelet & | operator= (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 ShapeletVector & | getValues () 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::BVec & | viewAsBVec () const |
View the shapelet as a BVec. More... | |
shapelet::BVec & | viewAsBVec () |
Private Attributes | |
ShapeletImpl * | pImpl |
Definition at line 50 of file Shapelet.h.
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.
typedef float lsst::meas::algorithms::Shapelet::PixelT |
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.
typedef Eigen::MatrixXd lsst::meas::algorithms::Shapelet::ShapeletCovariance |
Definition at line 98 of file Shapelet.h.
typedef Eigen::VectorXd lsst::meas::algorithms::Shapelet::ShapeletVector |
Definition at line 97 of file Shapelet.h.
Definition at line 100 of file Shapelet.h.
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.
order | the maximum value of p+q to measure. |
sigma | the scale size of the shapelet decomposition (arcsec). |
Definition at line 145 of file Shapelet.cc.
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.
order | the maximum value of p+q to measure. |
sigma | the scale size of the shapelet decomposition (arcsec). |
vector | the shapelet vector |
Definition at line 149 of file Shapelet.cc.
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.
order | the maximum value of p+q to measure. |
sigma | the scale size of the shapelet decomposition (arcsec). |
vector | the shapelet vector |
cov | the covariance matrix |
Definition at line 153 of file Shapelet.cc.
lsst::meas::algorithms::Shapelet::~Shapelet | ( | ) |
lsst::meas::algorithms::Shapelet::Shapelet | ( | const Shapelet & | rhs | ) |
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.
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.
double lsst::meas::algorithms::Shapelet::evaluateAt | ( | double | x, |
double | y | ||
) |
Definition at line 197 of file Shapelet.cc.
boost::shared_ptr< const Shapelet::ShapeletCovariance > lsst::meas::algorithms::Shapelet::getCovariance | ( | ) | const |
get the covariance matrix
Definition at line 185 of file Shapelet.cc.
int lsst::meas::algorithms::Shapelet::getOrder | ( | ) | const |
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.
double lsst::meas::algorithms::Shapelet::getSigma | ( | ) | const |
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.
bool lsst::meas::algorithms::Shapelet::hasCovariance | ( | ) | const |
does the shapelet have a covariance matrix stored?
Definition at line 182 of file Shapelet.cc.
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.
source | The source of which to measure the decomposition. |
pos | The position of the object |
isCentroidFixed | Is sigma fixed? or should it be allowed to vary? |
isSigmaFixed | Is sigma fixed? or should it be allowed to vary? |
aperture | The aperture size in arcsec for the measurement. |
exposure | The exposure on which to measure the decompoisition. |
okmask | The mask values that are ok to use |
Definition at line 358 of file Shapelet.cc.
void lsst::meas::algorithms::Shapelet::setSigma | ( | double | sigma | ) |
int lsst::meas::algorithms::Shapelet::size | ( | ) | const |
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.
shapelet::BVec & lsst::meas::algorithms::Shapelet::viewAsBVec | ( | ) |
Definition at line 411 of file Shapelet.cc.
|
private |
Definition at line 306 of file Shapelet.h.