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
Shapelet.h
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 #ifndef LSST_MEAS_ALGORITHMS_SHAPELET_H
3 #define LSST_MEAS_ALGORITHMS_SHAPELET_H
4 
5 /*
6  * LSST Data Management System
7  * Copyright 2008, 2009, 2010 LSST Corporation.
8  *
9  * This product includes software developed by the
10  * LSST Project (http://www.lsst.org/).
11  *
12  * This program is free software: you can redistribute it and/or modify
13  * it under the terms of the GNU General Public License as published by
14  * the Free Software Foundation, either version 3 of the License, or
15  * (at your option) any later version.
16  *
17  * This program is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  * GNU General Public License for more details.
21  *
22  * You should have received a copy of the LSST License Statement and
23  * the GNU General Public License along with this program. If not,
24  * see <http://www.lsstcorp.org/LegalNotices/>.
25  */
26 
34 #include <complex>
35 
36 #include "boost/shared_ptr.hpp"
37 #include "Eigen/Core"
38 
39 #include "lsst/afw/table/Source.h"
41 #include "lsst/afw/geom/Point.h"
42 
43 namespace lsst {
44 namespace meas {
45 namespace algorithms {
46 
47  class ShapeletImpl;
48  namespace shapelet { class BVec; }
49 
50  class Shapelet
51  {
92  public:
93  typedef float PixelT;
94  typedef boost::shared_ptr<Shapelet> Ptr;
95  typedef boost::shared_ptr<const Shapelet> ConstPtr;
96 
97  typedef Eigen::VectorXd ShapeletVector;
98  typedef Eigen::MatrixXd ShapeletCovariance;
99 
104 
116  Shapelet(
117  int order,
118  double sigma
119  );
120 
129  Shapelet(
130  int order,
131  double sigma,
132  const ShapeletVector& vector
133  );
134 
146  Shapelet(
147  int order,
148  double sigma,
149  const ShapeletVector& vector,
150  const ShapeletCovariance& cov
151  );
152 
156  ~Shapelet();
157 
161  Shapelet(const Shapelet& rhs);
162 
166  Shapelet& operator=(const Shapelet& rhs);
167 
171  int getOrder() const;
172 
176  double getSigma() const;
177 
181  int size() const;
182 
197  const ShapeletVector& getValues() const;
198 
202  bool hasCovariance() const;
203 
207  boost::shared_ptr<const ShapeletCovariance> getCovariance() const;
208 
212  void setSigma(double sigma);
213 
219  std::complex<double> getPQ(int p, int q);
220 
229  double evaluateAt(const PointD& pos);
230  double evaluateAt(double x, double y);
231 
275  bool measureFromImage(
276  const Source& source,
277  const PointD& pos,
278  bool isCentroidFixed,
279  bool isSigmaFixed,
280  double aperture,
281  const Exposure& exposure,
282  const MaskPixel okmask=0
283  );
284 
292  Shapelet(const shapelet::BVec& bvec);
293 
301  const shapelet::BVec& viewAsBVec() const;
303 
304  private :
305 
307  };
308 
325  Eigen::Matrix2d getJacobian(
326  const lsst::afw::image::Wcs& wcs,
327  const lsst::afw::geom::PointD& pos
328  );
329 
330 }}}
331 
332 #endif
int y
const ShapeletVector & getValues() const
get the values as a vector.
Definition: Shapelet.cc:179
A coordinate class intended to represent absolute positions.
~Shapelet()
Destructor needs to delete pImpl.
Definition: Shapelet.cc:161
void setSigma(double sigma)
set a new value of sigma
Definition: Shapelet.cc:188
A class to contain the data, WCS, and other information needed to describe an image of the sky...
Definition: Exposure.h:48
boost::uint16_t MaskPixel
double evaluateAt(const PointD &pos)
Evaluate f(x,y)
Definition: Shapelet.cc:194
tbl::Key< int > wcs
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
Definition: Shapelet.cc:358
lsst::afw::geom::PointD PointD
Definition: Shapelet.h:102
Implementation of the WCS standard for a any projection.
Definition: Wcs.h:107
double getSigma() const
get the scale size of the shapelet
Definition: Shapelet.cc:173
int getOrder() const
get the order of the shapelet
Definition: Shapelet.cc:170
const shapelet::BVec & viewAsBVec() const
View the shapelet as a BVec.
Definition: Shapelet.cc:410
Shapelet & operator=(const Shapelet &rhs)
op= does a deep copy
Definition: Shapelet.cc:164
boost::shared_ptr< const Shapelet > ConstPtr
Definition: Shapelet.h:95
Eigen::MatrixXd ShapeletCovariance
Definition: Shapelet.h:98
afw::table::Key< double > sigma
Definition: GaussianPsf.cc:43
lsst::afw::table::SourceRecord Source
Definition: Shapelet.h:100
bool hasCovariance() const
does the shapelet have a covariance matrix stored?
Definition: Shapelet.cc:182
boost::shared_ptr< const ShapeletCovariance > getCovariance() const
get the covariance matrix
Definition: Shapelet.cc:185
Eigen::Matrix2d getJacobian(const lsst::afw::image::Wcs &wcs, const lsst::afw::geom::PointD &pos)
a helper function to deal with the fact that Wcs doesn&#39;t directly return a Jacobian.
Definition: Shapelet.cc:209
boost::shared_ptr< Shapelet > Ptr
Definition: Shapelet.h:94
std::complex< double > getPQ(int p, int q)
Get a complex b_pq value.
Definition: Shapelet.cc:191
lsst::afw::image::Exposure< PixelT > Exposure
Definition: Shapelet.h:101
double x
lsst::afw::image::MaskPixel MaskPixel
Definition: Shapelet.h:103
Eigen::VectorXd ShapeletVector
Definition: Shapelet.h:97
Record class that contains measurements made on a single exposure.
Definition: Source.h:81
float PixelT
This class includes the basic functionality for shapelets.
Definition: Shapelet.h:93
int size() const
the size of the shapelet vector
Definition: Shapelet.cc:176
Shapelet(int order, double sigma)
Basic constructor requires order and sigma.
Definition: Shapelet.cc:145