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 | Static Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
lsst::afw::math::ChebyshevBoundedField Class Reference

A BoundedField based on 2-d Chebyshev polynomials of the first kind. More...

#include <ChebyshevBoundedField.h>

Inheritance diagram for lsst::afw::math::ChebyshevBoundedField:
lsst::afw::table::io::PersistableFacade< ChebyshevBoundedField > lsst::afw::math::BoundedField lsst::afw::table::io::PersistableFacade< BoundedField > lsst::afw::table::io::Persistable

Public Types

typedef
ChebyshevBoundedFieldControl 
Control
 

Public Member Functions

 ChebyshevBoundedField (afw::geom::Box2I const &bbox, ndarray::Array< double const, 2, 2 > const &coefficients)
 Initialize the field from its bounding box an coefficients. More...
 
ndarray::Array< double const, 2, 2 > getCoefficients () const
 Return the coefficient matrix. More...
 
boost::shared_ptr
< ChebyshevBoundedField
truncate (Control const &ctrl) const
 Return a new ChebyshevBoudedField with maximum orders set by the given control object. More...
 
boost::shared_ptr
< ChebyshevBoundedField
relocate (geom::Box2I const &bbox) const
 
virtual double evaluate (geom::Point2D const &position) const
 
virtual bool isPersistable () const
 ChebyshevBoundedField is always persistable. More...
 
virtual boost::shared_ptr
< BoundedField
operator* (double const scale) const
 
- Public Member Functions inherited from lsst::afw::math::BoundedField
double evaluate (double x, double y) const
 
ndarray::Array< double, 1, 1 > evaluate (ndarray::Array< double const, 1 > const &x, ndarray::Array< double const, 1 > const &y) const
 
geom::Box2I getBBox () const
 
template<typename T >
void fillImage (image::Image< T > &image, bool overlapOnly=false) const
 
template<typename T >
void addToImage (image::Image< T > &image, double scaleBy=1.0, bool overlapOnly=false) const
 
template<typename T >
void multiplyImage (image::Image< T > &image, bool overlapOnly=false) const
 
template<typename T >
void divideImage (image::Image< T > &image, bool overlapOnly=false) const
 
boost::shared_ptr< BoundedFieldoperator/ (double scale) const
 
virtual ~BoundedField ()
 
- Public Member Functions inherited from lsst::afw::table::io::Persistable
void writeFits (std::string const &fileName, std::string const &mode="w") const
 Write the object to a regular FITS file. More...
 
void writeFits (fits::MemFileManager &manager, std::string const &mode="w") const
 Write the object to a FITS image in memory. More...
 
void writeFits (fits::Fits &fitsfile) const
 Write the object to an already-open FITS object. More...
 
virtual ~Persistable ()
 

Static Public Member Functions

static boost::shared_ptr
< ChebyshevBoundedField
fit (afw::geom::Box2I const &bbox, ndarray::Array< double const, 1 > const &x, ndarray::Array< double const, 1 > const &y, ndarray::Array< double const, 1 > const &z, Control const &ctrl)
 Fit a Chebyshev approximation to non-gridded data with equal weights. More...
 
static boost::shared_ptr
< ChebyshevBoundedField
fit (afw::geom::Box2I const &bbox, ndarray::Array< double const, 1 > const &x, ndarray::Array< double const, 1 > const &y, ndarray::Array< double const, 1 > const &z, ndarray::Array< double const, 1 > const &w, Control const &ctrl)
 Fit a Chebyshev approximation to non-gridded data with unequal weights. More...
 
template<typename T >
static boost::shared_ptr
< ChebyshevBoundedField
fit (image::Image< T > const &image, Control const &ctrl)
 Fit a Chebyshev approximation to gridded data with equal weights. More...
 
- Static Public Member Functions inherited from lsst::afw::table::io::PersistableFacade< ChebyshevBoundedField >
static boost::shared_ptr
< ChebyshevBoundedField > 
readFits (fits::Fits &fitsfile)
 Read an object from an already open FITS object. More...
 
static boost::shared_ptr
< ChebyshevBoundedField > 
readFits (std::string const &fileName, int hdu=0)
 Read an object from a regular FITS file. More...
 
static boost::shared_ptr
< ChebyshevBoundedField > 
readFits (fits::MemFileManager &manager, int hdu=0)
 Read an object from a FITS file in memory. More...
 
- Static Public Member Functions inherited from lsst::afw::table::io::PersistableFacade< BoundedField >
static boost::shared_ptr
< BoundedField > 
readFits (fits::Fits &fitsfile)
 Read an object from an already open FITS object. More...
 
static boost::shared_ptr
< BoundedField > 
readFits (std::string const &fileName, int hdu=0)
 Read an object from a regular FITS file. More...
 
static boost::shared_ptr
< BoundedField > 
readFits (fits::MemFileManager &manager, int hdu=0)
 Read an object from a FITS file in memory. More...
 

Protected Member Functions

virtual std::string getPersistenceName () const
 Return the unique name used to persist this object and look up its factory. More...
 
virtual std::string getPythonModule () const
 Return the fully-qualified Python module that should be imported to guarantee that its factory is registered. More...
 
virtual void write (OutputArchiveHandle &handle) const
 Write the object to one or more catalogs. More...
 
- Protected Member Functions inherited from lsst::afw::math::BoundedField
 BoundedField (geom::Box2I const &bbox)
 
- Protected Member Functions inherited from lsst::afw::table::io::Persistable
 Persistable ()
 
 Persistable (Persistable const &other)
 
void operator= (Persistable const &other)
 

Private Member Functions

 ChebyshevBoundedField (afw::geom::Box2I const &bbox)
 

Private Attributes

geom::AffineTransform _toChebyshevRange
 
ndarray::Array< double const, 2, 2 > _coefficients
 

Additional Inherited Members

- Protected Types inherited from lsst::afw::table::io::Persistable
typedef io::OutputArchiveHandle OutputArchiveHandle
 

Detailed Description

A BoundedField based on 2-d Chebyshev polynomials of the first kind.

ChebyshevBoundedField supports fitting to gridded and non-gridded data, as well coefficient matrices with different x- and y-order.

There is currently quite a bit of duplication of functionality between ChebyshevBoundedField, ApproximateChebyshev, and Chebyshev1Function2; the intent is that ChebyshevBoundedField will ultimately replace ApproximateChebyshev and should be preferred over Chebyshev1Function2 when the parametrization interface that is part of the Function2 class is not needed.

Definition at line 69 of file ChebyshevBoundedField.h.

Member Typedef Documentation

Definition at line 75 of file ChebyshevBoundedField.h.

Constructor & Destructor Documentation

lsst::afw::math::ChebyshevBoundedField::ChebyshevBoundedField ( afw::geom::Box2I const &  bbox,
ndarray::Array< double const, 2, 2 > const &  coefficients 
)

Initialize the field from its bounding box an coefficients.

This constructor is mostly intended for testing purposes and persistence, but it also provides a way to initialize the object from Chebyshev coefficients derived from some external source.

Note that because the bounding box provided is always an integer bounding box, and LSST convention puts the center of each pixel at an integer, the actual floating-point domain of the Chebyshev functions is Box2D(bbox), that is, the box that contains the entirety of all the pixels included in the integer bounding box.

The coefficients are ordered [y,x], so the shape is (orderY+1, orderX+1), and the arguments to the Chebyshev functions are transformed such that the region Box2D(bbox) is mapped to [-1, 1]x[-1, 1].

Definition at line 58 of file ChebyshevBoundedField.cc.

61  : BoundedField(bbox),
62  _toChebyshevRange(makeChebyshevRangeTransform(geom::Box2D(bbox))),
63  _coefficients(coefficients)
64 {}
ndarray::Array< double const, 2, 2 > _coefficients
BoundedField(geom::Box2I const &bbox)
Definition: BoundedField.h:164
A floating-point coordinate rectangle geometry.
Definition: Box.h:271
lsst::afw::math::ChebyshevBoundedField::ChebyshevBoundedField ( afw::geom::Box2I const &  bbox)
explicitprivate

Definition at line 66 of file ChebyshevBoundedField.cc.

68  : BoundedField(bbox),
69  _toChebyshevRange(makeChebyshevRangeTransform(geom::Box2D(bbox)))
70 {}
BoundedField(geom::Box2I const &bbox)
Definition: BoundedField.h:164
A floating-point coordinate rectangle geometry.
Definition: Box.h:271

Member Function Documentation

double lsst::afw::math::ChebyshevBoundedField::evaluate ( geom::Point2D const &  position) const
virtual

Evaluate the field at the given point.

This is the only abstract method to be implemented by subclasses.

Subclasses should not provide bounds checking on the given position; this is the responsibility of the user, who can almost always do it more efficiently.

Implements lsst::afw::math::BoundedField.

Definition at line 310 of file ChebyshevBoundedField.cc.

310  {
311  geom::Point2D p = _toChebyshevRange(position);
312  return evaluateFunction1d(RecursionArrayImitator(_coefficients, p.getX()),
313  p.getY(), _coefficients.getSize<0>());
314 }
ndarray::Array< double const, 2, 2 > _coefficients
int getSize() const
Return the size of a specific dimension.
Definition: ArrayBase.h:126
boost::shared_ptr< ChebyshevBoundedField > lsst::afw::math::ChebyshevBoundedField::fit ( afw::geom::Box2I const &  bbox,
ndarray::Array< double const, 1 > const &  x,
ndarray::Array< double const, 1 > const &  y,
ndarray::Array< double const, 1 > const &  z,
Control const &  ctrl 
)
static

Fit a Chebyshev approximation to non-gridded data with equal weights.

Parameters
[in]bboxInteger bounding box of the resulting approximation. All given points must lie within Box2D(bbox).
[in]xArray of x coordinate values.
[in]yArray of y coordinate values.
[in]zArray of field values to be fit at each (x,y) point.
[in]ctrlSpecifies the orders and triangularity of the coefficient matrix.

Definition at line 161 of file ChebyshevBoundedField.cc.

167  {
168  // Initialize the result object, so we can make use of the AffineTransform it builds
170  // This packer object knows how to map the 2-d Chebyshev functions onto a 1-d array,
171  // using only those that the control says should have nonzero coefficients.
172  Packer const packer(ctrl);
173  // Create a "design matrix" for the linear least squares problem (A in min||Ax-b||)
174  ndarray::Array<double,2,2> matrix = makeMatrix(x, y, result->_toChebyshevRange, packer, ctrl);
175  // Solve the linear least squares problem.
176  LeastSquares lstsq = LeastSquares::fromDesignMatrix(matrix, z, LeastSquares::NORMAL_EIGENSYSTEM);
177  // Unpack the solution into a 2-d matrix, with zeros for values we didn't fit.
178  result->_coefficients = packer.unpack(lstsq.getSolution());
179  return result;
180 }
ndarray::Array< double const, 2, 2 > _coefficients
#define PTR(...)
Definition: base.h:41
ChebyshevBoundedField(afw::geom::Box2I const &bbox, ndarray::Array< double const, 2, 2 > const &coefficients)
Initialize the field from its bounding box an coefficients.
boost::shared_ptr< ChebyshevBoundedField > lsst::afw::math::ChebyshevBoundedField::fit ( afw::geom::Box2I const &  bbox,
ndarray::Array< double const, 1 > const &  x,
ndarray::Array< double const, 1 > const &  y,
ndarray::Array< double const, 1 > const &  z,
ndarray::Array< double const, 1 > const &  w,
Control const &  ctrl 
)
static

Fit a Chebyshev approximation to non-gridded data with unequal weights.

Parameters
[in]bboxInteger bounding box of the resulting approximation. All given points must lie within Box2D(bbox).
[in]xArray of x coordinate values.
[in]yArray of y coordinate values.
[in]zArray of field values to be fit at each (x,y) point.
[in]wArray of weights for each point in the fit. For points with Gaussian noise, w = 1/sigma.
[in]ctrlSpecifies the orders and triangularity of the coefficient matrix.

Definition at line 182 of file ChebyshevBoundedField.cc.

189  {
190  // Initialize the result object, so we can make use of the AffineTransform it builds
192  // This packer object knows how to map the 2-d Chebyshev functions onto a 1-d array,
193  // using only those that the control says should have nonzero coefficients.
194  Packer const packer(ctrl);
195  // Create a "design matrix" for the linear least squares problem ('A' in min||Ax-b||)
196  ndarray::Array<double,2,2> matrix = makeMatrix(x, y, result->_toChebyshevRange, packer, ctrl);
197  // We want to do weighted least squares, so we multiply both the data vector 'b' and the
198  // matrix 'A' by the weights.
199  matrix.asEigen<Eigen::ArrayXpr>().colwise() *= w.asEigen<Eigen::ArrayXpr>();
200  ndarray::Array<double,1,1> wz = ndarray::copy(z);
201  wz.asEigen<Eigen::ArrayXpr>() *= w.asEigen<Eigen::ArrayXpr>();
202  // Solve the linear least squares problem.
203  LeastSquares lstsq = LeastSquares::fromDesignMatrix(matrix, wz, LeastSquares::NORMAL_EIGENSYSTEM);
204  // Unpack the solution into a 2-d matrix, with zeros for values we didn't fit.
205  result->_coefficients = packer.unpack(lstsq.getSolution());
206  return result;
207 }
ndarray::Array< double const, 2, 2 > _coefficients
SelectEigenView< T >::Type copy(Eigen::EigenBase< T > const &other)
Copy an arbitrary Eigen expression into a new EigenView.
Definition: eigen.h:390
#define PTR(...)
Definition: base.h:41
ChebyshevBoundedField(afw::geom::Box2I const &bbox, ndarray::Array< double const, 2, 2 > const &coefficients)
Initialize the field from its bounding box an coefficients.
template<typename T >
boost::shared_ptr< ChebyshevBoundedField > lsst::afw::math::ChebyshevBoundedField::fit ( image::Image< T > const &  image,
Control const &  ctrl 
)
static

Fit a Chebyshev approximation to gridded data with equal weights.

Parameters
[in]imageThe Image containing the data to fit. image.getBBox(PARENT) is used as the bounding box of the BoundedField.
[in]ctrlSpecifies the orders and triangularity of the coefficient matrix.

Instantiated for float and double.

Note
if the image to be fit is a binned version of the actual image the field should correspond to, call relocate() with the unbinned image's bounding box after fitting.

Definition at line 210 of file ChebyshevBoundedField.cc.

213  {
214  // Initialize the result object, so we can make use of the AffineTransform it builds
215  geom::Box2I bbox = img.getBBox(image::PARENT);
217  // This packer object knows how to map the 2-d Chebyshev functions onto a 1-d array,
218  // using only those that the control says should have nonzero coefficients.
219  Packer const packer(ctrl);
220  ndarray::Array<double,2,2> matrix = makeMatrix(bbox, result->_toChebyshevRange, packer, ctrl);
221  // Flatten the data image into a 1-d vector.
222  ndarray::Array<double,2,2> imgCopy = ndarray::allocate(img.getArray().getShape());
223  imgCopy.deep() = img.getArray();
224  ndarray::Array<double const,1,1> z = ndarray::flatten<1>(imgCopy);
225  // Solve the linear least squares problem.
226  LeastSquares lstsq = LeastSquares::fromDesignMatrix(matrix, z, LeastSquares::NORMAL_EIGENSYSTEM);
227  // Unpack the solution into a 2-d matrix, with zeros for values we didn't fit.
228  result->_coefficients = packer.unpack(lstsq.getSolution());
229  return result;
230 }
ndarray::Array< double const, 2, 2 > _coefficients
#define PTR(...)
Definition: base.h:41
An integer coordinate rectangle.
Definition: Box.h:53
boost::enable_if_c< ((C+Nf-N)>=1), ArrayRef< T, Nf,(C+Nf-N)> >::type flatten(Array< T, N, C > const &input)
Create a view into an array with trailing contiguous dimensions merged.
Definition: casts.h:150
detail::SimpleInitializer< N > allocate(Vector< int, N > const &shape)
Create an expression that allocates uninitialized memory for an array.
ChebyshevBoundedField(afw::geom::Box2I const &bbox, ndarray::Array< double const, 2, 2 > const &coefficients)
Initialize the field from its bounding box an coefficients.
ndarray::Array<double const,2,2> lsst::afw::math::ChebyshevBoundedField::getCoefficients ( ) const
inline

Return the coefficient matrix.

The coefficients are ordered [y,x], so the shape is (orderY+1, orderX+1).

Definition at line 162 of file ChebyshevBoundedField.h.

162 { return _coefficients; }
ndarray::Array< double const, 2, 2 > _coefficients
std::string lsst::afw::math::ChebyshevBoundedField::getPersistenceName ( ) const
protectedvirtual

Return the unique name used to persist this object and look up its factory.

Must be less than ArchiveIndexSchema::MAX_NAME_LENGTH characters.

Reimplemented from lsst::afw::table::io::Persistable.

Definition at line 379 of file ChebyshevBoundedField.cc.

379  {
380  return getChebyshevBoundedFieldPersistenceName();
381 }
std::string lsst::afw::math::ChebyshevBoundedField::getPythonModule ( ) const
protectedvirtual

Return the fully-qualified Python module that should be imported to guarantee that its factory is registered.

Must be less than ArchiveIndexSchema::MAX_MODULE_LENGTH characters.

Will be ignored if empty.

Reimplemented from lsst::afw::table::io::Persistable.

Definition at line 383 of file ChebyshevBoundedField.cc.

383  {
384  return "lsst.afw.math";
385 }
virtual bool lsst::afw::math::ChebyshevBoundedField::isPersistable ( ) const
inlinevirtual

ChebyshevBoundedField is always persistable.

Reimplemented from lsst::afw::table::io::Persistable.

Definition at line 181 of file ChebyshevBoundedField.h.

181 { return true; }
boost::shared_ptr< BoundedField > lsst::afw::math::ChebyshevBoundedField::operator* ( double const  scale) const
virtual

Return a scaled BoundedField

Parameters
[in]scaleScaling factor

Implements lsst::afw::math::BoundedField.

Definition at line 398 of file ChebyshevBoundedField.cc.

398  {
399  return boost::make_shared<ChebyshevBoundedField>(getBBox(), ndarray::copy(getCoefficients()*scale));
400 }
SelectEigenView< T >::Type copy(Eigen::EigenBase< T > const &other)
Copy an arbitrary Eigen expression into a new EigenView.
Definition: eigen.h:390
geom::Box2I getBBox() const
Definition: BoundedField.h:95
ndarray::Array< double const, 2, 2 > getCoefficients() const
Return the coefficient matrix.
boost::shared_ptr< ChebyshevBoundedField > lsst::afw::math::ChebyshevBoundedField::relocate ( geom::Box2I const &  bbox) const

Return a new ChebyshevBoundedField with domain set to the given bounding box.

Because this leaves the coefficients unchanged, it is equivalent to transforming the function by the affine transform that maps the old box to the new one.

Definition at line 261 of file ChebyshevBoundedField.cc.

261  {
262  return boost::make_shared<ChebyshevBoundedField>(bbox, _coefficients);
263 }
ndarray::Array< double const, 2, 2 > _coefficients
boost::shared_ptr< ChebyshevBoundedField > lsst::afw::math::ChebyshevBoundedField::truncate ( Control const &  ctrl) const

Return a new ChebyshevBoudedField with maximum orders set by the given control object.

Definition at line 234 of file ChebyshevBoundedField.cc.

234  {
235  if (ctrl.orderX >= _coefficients.getSize<1>()) {
236  throw LSST_EXCEPT(
237  pex::exceptions::LengthError,
238  (boost::format("New x order (%d) exceeds old x order (%d)")
239  % ctrl.orderX % (_coefficients.getSize<1>() - 1)).str()
240  );
241  }
242  if (ctrl.orderY >= _coefficients.getSize<0>()) {
243  throw LSST_EXCEPT(
244  pex::exceptions::LengthError,
245  (boost::format("New y order (%d) exceeds old y order (%d)")
246  % ctrl.orderY % (_coefficients.getSize<0>() - 1)).str()
247  );
248  }
249  ndarray::Array<double,2,2> coefficients = ndarray::allocate(ctrl.orderY + 1, ctrl.orderX + 1);
250  coefficients.deep()
251  = _coefficients[ndarray::view(0, ctrl.orderY + 1)(0, ctrl.orderX + 1)];
252  if (ctrl.triangular) {
253  Packer packer(ctrl);
254  ndarray::Array<double,1,1> packed = ndarray::allocate(packer.size);
255  packer.pack(packed, coefficients);
256  packer.unpack(coefficients, packed);
257  }
258  return boost::make_shared<ChebyshevBoundedField>(getBBox(), coefficients);
259 }
View< boost::fusion::vector1< index::Full > > view()
Start a view definition that includes the entire first dimension.
Definition: views.h:131
ndarray::Array< double const, 2, 2 > _coefficients
int getSize() const
Return the size of a specific dimension.
Definition: ArrayBase.h:126
geom::Box2I getBBox() const
Definition: BoundedField.h:95
detail::SimpleInitializer< N > allocate(Vector< int, N > const &shape)
Create an expression that allocates uninitialized memory for an array.
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46
Deep const deep() const
Return an ArrayRef view to this.
Definition: ArrayBase.h:178
ndarray::Array< double const, 2, 2 > coefficients
void lsst::afw::math::ChebyshevBoundedField::write ( OutputArchiveHandle handle) const
protectedvirtual

Write the object to one or more catalogs.

The handle object passed to this function provides an interface for adding new catalogs and adding nested objects to the same archive (while checking for duplicates). See OutputArchiveHandle for more information.

Reimplemented from lsst::afw::table::io::Persistable.

Definition at line 387 of file ChebyshevBoundedField.cc.

387  {
388  PersistenceHelper const keys(_coefficients.getSize<1>(), _coefficients.getSize<0>());
389  table::BaseCatalog catalog = handle.makeCatalog(keys.schema);
390  PTR(table::BaseRecord) record = catalog.addNew();
391  record->set(keys.orderX, _coefficients.getSize<1>() - 1);
392  record->set(keys.bboxMin, getBBox().getMin());
393  record->set(keys.bboxMax, getBBox().getMax());
394  (*record)[keys.coefficients].deep() = ndarray::flatten<1>(_coefficients);
395  handle.saveCatalog(catalog);
396 }
table::PointKey< int > bboxMax
ndarray::Array< double const, 2, 2 > _coefficients
CatalogT< BaseRecord > BaseCatalog
Definition: fwd.h:61
#define PTR(...)
Definition: base.h:41
int getSize() const
Return the size of a specific dimension.
Definition: ArrayBase.h:126
geom::Box2I getBBox() const
Definition: BoundedField.h:95
boost::enable_if_c< ((C+Nf-N)>=1), ArrayRef< T, Nf,(C+Nf-N)> >::type flatten(Array< T, N, C > const &input)
Create a view into an array with trailing contiguous dimensions merged.
Definition: casts.h:150
table::PointKey< int > bboxMin
table::Key< int > orderX
ndarray::Array< double const, 2, 2 > coefficients

Member Data Documentation

ndarray::Array<double const,2,2> lsst::afw::math::ChebyshevBoundedField::_coefficients
private

Definition at line 201 of file ChebyshevBoundedField.h.

geom::AffineTransform lsst::afw::math::ChebyshevBoundedField::_toChebyshevRange
private

Definition at line 200 of file ChebyshevBoundedField.h.


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