LSST Applications g0000d66e7c+ce78115f25,g0485b4d2cb+c8d56b10d4,g0fba68d861+fcbc158cd0,g1ec0fe41b4+3e153770da,g1fd858c14a+57ee4e1624,g2440f9efcc+8c5ae1fdc5,g35bb328faa+8c5ae1fdc5,g4d2262a081+1e04cc5a47,g53246c7159+8c5ae1fdc5,g55585698de+7a33f081c8,g56a49b3a55+b9d5cac73f,g60b5630c4e+7a33f081c8,g67b6fd64d1+035c836e50,g78460c75b0+7e33a9eb6d,g786e29fd12+668abc6043,g7ac00fbb6c+b938379438,g8352419a5c+8c5ae1fdc5,g8852436030+5ba78a36c9,g89139ef638+035c836e50,g94187f82dc+7a33f081c8,g989de1cb63+035c836e50,g9d31334357+7a33f081c8,g9f33ca652e+e34120223a,ga815be3f0b+911242149a,gabe3b4be73+8856018cbb,gabf8522325+21619da9f3,gb1101e3267+0b44b44611,gb89ab40317+035c836e50,gc91f06edcd+e59fb3c9bc,gcf25f946ba+5ba78a36c9,gd6cbbdb0b4+958adf5c1f,gde0f65d7ad+6c98dcc924,ge278dab8ac+83c63f4893,ge410e46f29+035c836e50,gf35d7ec915+97dd712d81,gf5e32f922b+8c5ae1fdc5,gf67bdafdda+035c836e50,gf6800124b1+1714c04baa,w.2025.19
LSST Data Management Base Package
Loading...
Searching...
No Matches
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

using Control = ChebyshevBoundedFieldControl
 

Public Member Functions

 ChebyshevBoundedField (lsst::geom::Box2I const &bbox, ndarray::Array< double const, 2, 2 > const &coefficients)
 Initialize the field from its bounding box an coefficients.
 
 ChebyshevBoundedField (ChebyshevBoundedField const &)
 
 ChebyshevBoundedField (ChebyshevBoundedField &&)
 
ChebyshevBoundedFieldoperator= (ChebyshevBoundedField const &)=delete
 
ChebyshevBoundedFieldoperator= (ChebyshevBoundedField &&)=delete
 
 ~ChebyshevBoundedField () override
 
ndarray::Array< double const, 2, 2 > getCoefficients () const
 Return the coefficient matrix.
 
std::shared_ptr< ChebyshevBoundedFieldtruncate (Control const &ctrl) const
 Return a new ChebyshevBoudedField with maximum orders set by the given control object.
 
std::shared_ptr< ChebyshevBoundedFieldrelocate (lsst::geom::Box2I const &bbox) const
 Return a new ChebyshevBoundedField with domain set to the given bounding box.
 
double evaluate (lsst::geom::Point2D const &position) const override
 Evaluate the field at the given point.
 
double integrate () const override
 Compute the integral of this function over its bounding-box.
 
double mean () const override
 Compute the mean of this function over its bounding-box.
 
bool isPersistable () const noexcept override
 ChebyshevBoundedField is always persistable.
 
std::shared_ptr< BoundedFieldoperator* (double const scale) const override
 Return a scaled BoundedField.
 
bool operator== (BoundedField const &rhs) const override
 BoundedFields (of the same sublcass) are equal if their bounding boxes and parameters are equal.
 
double evaluate (double x, double y) const
 Evaluate the field at the given point.
 
virtual ndarray::Array< double, 1, 1 > evaluate (ndarray::Array< double const, 1 > const &x, ndarray::Array< double const, 1 > const &y) const
 Evaluate the field at multiple arbitrary points.
 
lsst::geom::Box2I getBBox () const
 Return the bounding box that defines the region where the field is valid.
 
template<typename T>
void fillImage (image::Image< T > &image, bool overlapOnly=false, int xStep=1, int yStep=1) const
 Assign the field to an image, overwriting values already present.
 
template<typename T>
void addToImage (image::Image< T > &image, double scaleBy=1.0, bool overlapOnly=false, int xStep=1, int yStep=1) const
 Add the field or a constant multiple of it to an image in-place.
 
template<typename T>
void multiplyImage (image::Image< T > &image, bool overlapOnly=false, int xStep=1, int yStep=1) const
 Multiply an image by the field in-place.
 
template<typename T>
void multiplyImage (image::MaskedImage< T > &image, bool overlapOnly=false) const
 
template<typename T>
void divideImage (image::Image< T > &image, bool overlapOnly=false, int xStep=1, int yStep=1) const
 Divide an image by the field in-place.
 
template<typename T>
void divideImage (image::MaskedImage< T > &image, bool overlapOnly=false) const
 
std::shared_ptr< BoundedFieldoperator/ (double scale) const
 
bool operator!= (BoundedField const &rhs) const
 BoundedFields (of the same sublcass) are equal if their bounding boxes and parameters are equal.
 
void writeFits (std::string const &fileName, std::string const &mode="w") const
 Write the object to a regular FITS file.
 
void writeFits (fits::MemFileManager &manager, std::string const &mode="w") const
 Write the object to a FITS image in memory.
 
void writeFits (fits::Fits &fitsfile) const
 Write the object to an already-open FITS object.
 

Static Public Member Functions

static std::shared_ptr< ChebyshevBoundedFieldfit (lsst::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.
 
static std::shared_ptr< ChebyshevBoundedFieldfit (lsst::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.
 
static ndarray::Array< double, 2, 2 > makeFitMatrix (lsst::geom::Box2I const &bbox, ndarray::Array< double const, 1 > const &x, ndarray::Array< double const, 1 > const &y, Control const &ctrl)
 Construct a matrix that can be used to fit or evaluate a Chebyshev polynomial at a set of points.
 
template<typename T>
static std::shared_ptr< ChebyshevBoundedFieldfit (image::Image< T > const &image, Control const &ctrl)
 Fit a Chebyshev approximation to gridded data with equal weights.
 
static std::shared_ptr< ChebyshevBoundedFieldreadFits (fits::Fits &fitsfile)
 Read an object from an already open FITS object.
 
static std::shared_ptr< ChebyshevBoundedFieldreadFits (std::string const &fileName, int hdu=fits::DEFAULT_HDU)
 Read an object from a regular FITS file.
 
static std::shared_ptr< ChebyshevBoundedFieldreadFits (fits::MemFileManager &manager, int hdu=fits::DEFAULT_HDU)
 Read an object from a FITS file in memory.
 
static std::shared_ptr< ChebyshevBoundedFielddynamicCast (std::shared_ptr< Persistable > const &ptr)
 Dynamically cast a shared_ptr.
 
static std::shared_ptr< BoundedFieldreadFits (fits::Fits &fitsfile)
 Read an object from an already open FITS object.
 
static std::shared_ptr< BoundedFieldreadFits (std::string const &fileName, int hdu=fits::DEFAULT_HDU)
 Read an object from a regular FITS file.
 
static std::shared_ptr< BoundedFieldreadFits (fits::MemFileManager &manager, int hdu=fits::DEFAULT_HDU)
 Read an object from a FITS file in memory.
 
static std::shared_ptr< BoundedFielddynamicCast (std::shared_ptr< Persistable > const &ptr)
 Dynamically cast a shared_ptr.
 

Protected Types

using OutputArchiveHandle = io::OutputArchiveHandle
 

Protected Member Functions

std::string getPersistenceName () const override
 Return the unique name used to persist this object and look up its factory.
 
std::string getPythonModule () const override
 Return the fully-qualified Python module that should be imported to guarantee that its factory is registered.
 
void write (OutputArchiveHandle &handle) const override
 Write the object to one or more catalogs.
 

Private Member Functions

std::string toString () const override
 

Detailed Description

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

The 2-d Chebyshev polynomial used here is defined as:

\[f(x,y) = \sum_i \sum_j a_{i,j} T_i(x) T_j(y) \]

where \(T_n(x)\) is the n-th order Chebyshev polynomial of \(x\) and \(a_{i,j}\) is the corresponding coefficient of the (i,j) polynomial term.

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 76 of file ChebyshevBoundedField.h.

Member Typedef Documentation

◆ Control

◆ OutputArchiveHandle

using lsst::afw::table::io::Persistable::OutputArchiveHandle = io::OutputArchiveHandle
protectedinherited

Definition at line 108 of file Persistable.h.

Constructor & Destructor Documentation

◆ ChebyshevBoundedField() [1/3]

lsst::afw::math::ChebyshevBoundedField::ChebyshevBoundedField ( lsst::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 lsst::geom::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 lsst::geom::Box2D(bbox) is mapped to [-1, 1]x[-1, 1].

Example:

bbox = lsst::geom::Box2I(lsst::geom::Point2I(10, 20), lsst::geom::Point2I(30, 40));
ndarray::Array<double, 2, 2> coeffs = ndarray::allocate(ndarray::makeVector(2, 2));
coeffs[0][0] = 1;
coeffs[1][0] = 2;
coeffs[0][1] = 3;
coeffs[1][1] = 4;
ndarray::Array<double, 2, 2> coeffs = ndarray::external(data);
poly = ChebyshevBoundedField(bbox, coeffs);

will result in the following polynomial:

\[f(x,y) = 1 T_0(x) T_0(y) + 2 T_0(x) T_1(y) + 3 T_1(x) T_0(y) + 4 T_1(x) T_1(y) \]

Definition at line 59 of file ChebyshevBoundedField.cc.

61 : BoundedField(bbox),
62 _toChebyshevRange(makeChebyshevRangeTransform(lsst::geom::Box2D(bbox))),
63 _coefficients(coefficients) {}
BoundedField(BoundedField const &)=default

◆ ChebyshevBoundedField() [2/3]

lsst::afw::math::ChebyshevBoundedField::ChebyshevBoundedField ( ChebyshevBoundedField const & )
default

◆ ChebyshevBoundedField() [3/3]

lsst::afw::math::ChebyshevBoundedField::ChebyshevBoundedField ( ChebyshevBoundedField && )
default

◆ ~ChebyshevBoundedField()

lsst::afw::math::ChebyshevBoundedField::~ChebyshevBoundedField ( )
overridedefault

Member Function Documentation

◆ addToImage()

template<typename T>
template void lsst::afw::math::BoundedField::addToImage ( image::Image< T > & image,
double scaleBy = 1.0,
bool overlapOnly = false,
int xStep = 1,
int yStep = 1 ) const
inherited

Add the field or a constant multiple of it to an image in-place.

Parameters
[out]imageImage to add to.
[in]scaleByMultiply the field by this before adding it to the image.
[in]overlapOnlyIf true, only modify the region in the intersection of image.getBBox(image::PARENT) and this->getBBox().
[in]xStepDistance between grid points in X to evaluate; values between grid points will be linearly interpolated.
[in]yStepDistance between grid points in Y to evaluate; values between grid points will be linearly interpolated.
Exceptions
pex::exceptions::RuntimeErrorif the bounding boxes do not overlap and overlapOnly=false.

Definition at line 329 of file BoundedField.cc.

330 {
331 applyToImage(*this, img, ScaledAdd(scaleBy), overlapOnly, xStep, yStep);
332}

◆ divideImage() [1/2]

template<typename T>
void lsst::afw::math::BoundedField::divideImage ( image::Image< T > & image,
bool overlapOnly = false,
int xStep = 1,
int yStep = 1 ) const
inherited

Divide an image by the field in-place.

Parameters
[out]imageImage to fill.
[in]overlapOnlyIf true, only modify the region in the intersection of image.getBBox(image::PARENT) and this->getBBox().
[in]xStepDistance between grid points in X to evaluate; values between grid points will be linearly interpolated.
[in]yStepDistance between grid points in Y to evaluate; values between grid points will be linearly interpolated.
Exceptions
pex::exceptions::RuntimeErrorif the bounding boxes do not overlap and overlapOnly=false.

Definition at line 345 of file BoundedField.cc.

345 {
346 applyToImage(*this, img, Divide(), overlapOnly, xStep, yStep);
347}

◆ divideImage() [2/2]

template<typename T>
void lsst::afw::math::BoundedField::divideImage ( image::MaskedImage< T > & image,
bool overlapOnly = false ) const
inherited

Definition at line 350 of file BoundedField.cc.

350 {
351 applyToImage(*this, img, Divide(), overlapOnly, 1, 1);
352}

◆ dynamicCast() [1/2]

Dynamically cast a shared_ptr.

Dynamically cast a shared pointer and raise on failure.

You must provide an explicit template instantiation in the .cc file for each class that inherits from PersistableFacade. Designed to work around RTTI issues on macOS with hidden symbols;

Exceptions
lsst::pex::exceptions::LogicErrorif the cast fails

param[in] ptr The pointer to be cast.

Returns
The cast pointer.
Exceptions
lsst::pex::exceptions::TypeErrorIf the dynamic cast fails.

Definition at line 218 of file Persistable.cc.

18 {
20 if (!result) {
21 throw LSST_EXCEPT(pex::exceptions::TypeError, "Dynamic pointer cast failed");
22 }
23 return result;
24}
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition Exception.h:48
A CRTP facade class for subclasses of Persistable.
T dynamic_pointer_cast(T... args)

◆ dynamicCast() [2/2]

Dynamically cast a shared_ptr.

Dynamically cast a shared pointer and raise on failure.

You must provide an explicit template instantiation in the .cc file for each class that inherits from PersistableFacade. Designed to work around RTTI issues on macOS with hidden symbols;

Exceptions
lsst::pex::exceptions::LogicErrorif the cast fails

param[in] ptr The pointer to be cast.

Returns
The cast pointer.
Exceptions
lsst::pex::exceptions::TypeErrorIf the dynamic cast fails.

Definition at line 218 of file Persistable.cc.

18 {
20 if (!result) {
21 throw LSST_EXCEPT(pex::exceptions::TypeError, "Dynamic pointer cast failed");
22 }
23 return result;
24}

◆ evaluate() [1/3]

double lsst::afw::math::BoundedField::evaluate ( double x,
double y ) const
inline

Evaluate the field at the given point.

This delegates to the evaluate() method that takes lsst::geom::Point2D.

There is no bounds-checking on the given position; this is the responsibility of the user, who can almost always do it more efficiently.

Definition at line 76 of file BoundedField.h.

76{ return evaluate(lsst::geom::Point2D(x, y)); }
double evaluate(lsst::geom::Point2D const &position) const override
Evaluate the field at the given point.
Point< double, 2 > Point2D
Definition Point.h:324

◆ evaluate() [2/3]

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

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 296 of file ChebyshevBoundedField.cc.

296 {
297 lsst::geom::Point2D p = _toChebyshevRange(position);
298 return evaluateFunction1d(RecursionArrayImitator(_coefficients, p.getX()), p.getY(),
299 _coefficients.getSize<0>());
300}

◆ evaluate() [3/3]

ndarray::Array< double, 1, 1 > lsst::afw::math::BoundedField::evaluate ( ndarray::Array< double const, 1 > const & x,
ndarray::Array< double const, 1 > const & y ) const
virtual

Evaluate the field at multiple arbitrary points.

Parameters
[in]xarray of x coordinates, same shape as y
[in]yarray of y coordinates, same shape as x
Returns
an array of output values, same shape as x and y

There is no bounds-checking on the given positions; this is the responsibility of the user, who can almost always do it more efficiently.

Reimplemented from lsst::afw::math::BoundedField.

Definition at line 88 of file BoundedField.cc.

41 {
42 ndarray::Array<double, 1, 1> out = ndarray::allocate(x.getSize<0>());
43 for (int i = 0, n = x.getSize<0>(); i < n; ++i) {
44 out[i] = evaluate(x[i], y[i]);
45 }
46 return out;
47}

◆ fillImage()

template<typename T>
template void lsst::afw::math::BoundedField::fillImage ( image::Image< T > & image,
bool overlapOnly = false,
int xStep = 1,
int yStep = 1 ) const
inherited

Assign the field to an image, overwriting values already present.

Parameters
[out]imageImage to fill.
[in]overlapOnlyIf true, only modify the region in the intersection of image.getBBox(image::PARENT) and this->getBBox().
[in]xStepDistance between grid points in X to evaluate; values between grid points will be linearly interpolated.
[in]yStepDistance between grid points in Y to evaluate; values between grid points will be linearly interpolated.
Exceptions
pex::exceptions::RuntimeErrorif the bounding boxes do not overlap and overlapOnly=false.

Definition at line 324 of file BoundedField.cc.

324 {
325 applyToImage(*this, img, Assign(), overlapOnly, xStep, yStep);
326}

◆ fit() [1/3]

template<typename T>
std::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 208 of file ChebyshevBoundedField.cc.

209 {
210 // Initialize the result object, so we can make use of the AffineTransform it builds
211 lsst::geom::Box2I bbox = img.getBBox(image::PARENT);
212 std::shared_ptr<ChebyshevBoundedField> result(new ChebyshevBoundedField(bbox));
213 // This packer object knows how to map the 2-d Chebyshev functions onto a 1-d array,
214 // using only those that the control says should have nonzero coefficients.
215 Packer const packer(ctrl);
216 ndarray::Array<double, 2, 2> matrix = makeMatrix(bbox, result->_toChebyshevRange, packer, ctrl);
217 // Flatten the data image into a 1-d vector.
218 ndarray::Array<double, 2, 2> imgCopy = ndarray::allocate(img.getArray().getShape());
219 imgCopy.deep() = img.getArray();
220 ndarray::Array<double const, 1, 1> z = ndarray::flatten<1>(imgCopy);
221 // Solve the linear least squares problem.
223 // Unpack the solution into a 2-d matrix, with zeros for values we didn't fit.
224 result->_coefficients = packer.unpack(lstsq.getSolution());
225 return result;
226}
ChebyshevBoundedField(lsst::geom::Box2I const &bbox, ndarray::Array< double const, 2, 2 > const &coefficients)
Initialize the field from its bounding box an coefficients.
static LeastSquares fromDesignMatrix(ndarray::Array< T1, 2, C1 > const &design, ndarray::Array< T2, 1, C2 > const &data, Factorization factorization=NORMAL_EIGENSYSTEM)
Initialize from the design matrix and data vector given as ndarrays.
@ NORMAL_EIGENSYSTEM
Use the normal equations with a symmetric Eigensystem decomposition.

◆ fit() [2/3]

std::shared_ptr< ChebyshevBoundedField > lsst::afw::math::ChebyshevBoundedField::fit ( lsst::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 lsst::geom::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 148 of file ChebyshevBoundedField.cc.

152 {
153 // Initialize the result object, so we can make use of the AffineTransform it builds
154 std::shared_ptr<ChebyshevBoundedField> result(new ChebyshevBoundedField(bbox));
155 // This packer object knows how to map the 2-d Chebyshev functions onto a 1-d array,
156 // using only those that the control says should have nonzero coefficients.
157 Packer const packer(ctrl);
158 // Create a "design matrix" for the linear least squares problem (A in min||Ax-b||)
159 ndarray::Array<double, 2, 2> matrix = makeMatrix(x, y, result->_toChebyshevRange, packer, ctrl);
160 // Solve the linear least squares problem.
162 // Unpack the solution into a 2-d matrix, with zeros for values we didn't fit.
163 result->_coefficients = packer.unpack(lstsq.getSolution());
164 return result;
165}

◆ fit() [3/3]

std::shared_ptr< ChebyshevBoundedField > lsst::afw::math::ChebyshevBoundedField::fit ( lsst::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 lsst::geom::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 167 of file ChebyshevBoundedField.cc.

172 {
173 // Initialize the result object, so we can make use of the AffineTransform it builds
174 std::shared_ptr<ChebyshevBoundedField> result(new ChebyshevBoundedField(bbox));
175 // This packer object knows how to map the 2-d Chebyshev functions onto a 1-d array,
176 // using only those that the control says should have nonzero coefficients.
177 Packer const packer(ctrl);
178 // Create a "design matrix" for the linear least squares problem ('A' in min||Ax-b||)
179 ndarray::Array<double, 2, 2> matrix = makeMatrix(x, y, result->_toChebyshevRange, packer, ctrl);
180 // We want to do weighted least squares, so we multiply both the data vector 'b' and the
181 // matrix 'A' by the weights.
182 ndarray::asEigenArray(matrix).colwise() *= ndarray::asEigenArray(w);
183 ndarray::Array<double, 1, 1> wz = ndarray::copy(z);
184 ndarray::asEigenArray(wz) *= ndarray::asEigenArray(w);
185 // Solve the linear least squares problem.
187 // Unpack the solution into a 2-d matrix, with zeros for values we didn't fit.
188 result->_coefficients = packer.unpack(lstsq.getSolution());
189 return result;
190}

◆ getBBox()

lsst::geom::Box2I lsst::afw::math::BoundedField::getBBox ( ) const
inlineinherited

Return the bounding box that defines the region where the field is valid.

Because this is an integer bounding box, its minimum and maximum positions are the centers of the pixels where the field is valid, but the field can be assumed to be valid to the edges of those pixels, which is the boundary you'd get by converting the returned lsst::geom::Box2I into a lsst::geom::Box2D.

Definition at line 113 of file BoundedField.h.

113{ return _bbox; }

◆ getCoefficients()

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 203 of file ChebyshevBoundedField.h.

203{ return _coefficients; }

◆ getPersistenceName()

std::string lsst::afw::math::ChebyshevBoundedField::getPersistenceName ( ) const
overrideprotectedvirtual

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 372 of file ChebyshevBoundedField.cc.

372 {
373 return getChebyshevBoundedFieldPersistenceName();
374}

◆ getPythonModule()

std::string lsst::afw::math::ChebyshevBoundedField::getPythonModule ( ) const
overrideprotectedvirtual

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 376 of file ChebyshevBoundedField.cc.

376{ return "lsst.afw.math"; }

◆ integrate()

double lsst::afw::math::ChebyshevBoundedField::integrate ( ) const
overridevirtual

Compute the integral of this function over its bounding-box.

Returns
The value of the integral.

Reimplemented from lsst::afw::math::BoundedField.

Definition at line 311 of file ChebyshevBoundedField.cc.

311 {
312 double result = 0;
313 double determinant = getBBox().getArea() / 4.0;
314 for (ndarray::Size j = 0; j < _coefficients.getSize<0>(); j++) {
315 for (ndarray::Size i = 0; i < _coefficients.getSize<1>(); i++) {
316 result += _coefficients[j][i] * integrateTn(i) * integrateTn(j);
317 }
318 }
319 return result * determinant;
320}
lsst::geom::Box2I getBBox() const
Return the bounding box that defines the region where the field is valid.
int getArea() const
Definition Box.h:189

◆ isPersistable()

bool lsst::afw::math::ChebyshevBoundedField::isPersistable ( ) const
inlineoverridevirtualnoexcept

ChebyshevBoundedField is always persistable.

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

Definition at line 228 of file ChebyshevBoundedField.h.

228{ return true; }

◆ makeFitMatrix()

ndarray::Array< double, 2, 2 > lsst::afw::math::ChebyshevBoundedField::makeFitMatrix ( lsst::geom::Box2I const & bbox,
ndarray::Array< double const, 1 > const & x,
ndarray::Array< double const, 1 > const & y,
Control const & ctrl )
static

Construct a matrix that can be used to fit or evaluate a Chebyshev polynomial at a set of points.

Parameters
[in]bboxInteger bounding box of the resulting approximation. All given points must lie within lsst::geom::Box2D(bbox).
[in]xArray of x coordinate values.
[in]yArray of y coordinate values (same size as x).
[in]ctrlSpecifies the orders and triangularity of the coefficient matrix.
Returns
A matrix with shape (x.size(), ctrl.computeSize()). The order of coefficients is unspecified, but will always be the same for the same ctrl (so you can make one matrix to fit coefficients from some points, and another to evaluate the polynomial at other points).

Definition at line 192 of file ChebyshevBoundedField.cc.

197 {
198 // Initialize a temporary result object, so we can make use of the AffineTransform it builds
199 std::shared_ptr<ChebyshevBoundedField> result(new ChebyshevBoundedField(bbox));
200 // This packer object knows how to map the 2-d Chebyshev functions onto a 1-d array,
201 // using only those that the control says should have nonzero coefficients.
202 Packer const packer(ctrl);
203 // Create a "design matrix" for the linear least squares problem (A in min||Ax-b||)
204 return makeMatrix(x, y, result->_toChebyshevRange, packer, ctrl);
205}

◆ mean()

double lsst::afw::math::ChebyshevBoundedField::mean ( ) const
overridevirtual

Compute the mean of this function over its bounding-box.

Returns
The value of the mean.

Reimplemented from lsst::afw::math::BoundedField.

Definition at line 322 of file ChebyshevBoundedField.cc.

322{ return integrate() / getBBox().getArea(); }
double integrate() const override
Compute the integral of this function over its bounding-box.

◆ multiplyImage() [1/2]

template<typename T>
void lsst::afw::math::BoundedField::multiplyImage ( image::Image< T > & image,
bool overlapOnly = false,
int xStep = 1,
int yStep = 1 ) const
inherited

Multiply an image by the field in-place.

Parameters
[out]imageImage to fill.
[in]overlapOnlyIf true, only modify the region in the intersection of image.getBBox(image::PARENT) and this->getBBox().
[in]xStepDistance between grid points in X to evaluate; values between grid points will be linearly interpolated.
[in]yStepDistance between grid points in Y to evaluate; values between grid points will be linearly interpolated.
Exceptions
pex::exceptions::RuntimeErrorif the bounding boxes do not overlap and overlapOnly=false.

Definition at line 335 of file BoundedField.cc.

335 {
336 applyToImage(*this, img, Multiply(), overlapOnly, xStep, yStep);
337}

◆ multiplyImage() [2/2]

template<typename T>
void lsst::afw::math::BoundedField::multiplyImage ( image::MaskedImage< T > & image,
bool overlapOnly = false ) const
inherited

Definition at line 340 of file BoundedField.cc.

340 {
341 applyToImage(*this, img, Multiply(), overlapOnly, 1, 1);
342}

◆ operator!=()

bool lsst::afw::math::BoundedField::operator!= ( BoundedField const & rhs) const
inlineinherited

BoundedFields (of the same sublcass) are equal if their bounding boxes and parameters are equal.

Definition at line 206 of file BoundedField.h.

206{ return !(*this == rhs); };

◆ operator*()

std::shared_ptr< BoundedField > lsst::afw::math::ChebyshevBoundedField::operator* ( double const scale) const
overridevirtual

Return a scaled BoundedField.

Parameters
[in]scaleScaling factor

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

Definition at line 390 of file ChebyshevBoundedField.cc.

390 {
391 return std::make_shared<ChebyshevBoundedField>(getBBox(), ndarray::copy(getCoefficients() * scale));
392}
ndarray::Array< double const, 2, 2 > getCoefficients() const
Return the coefficient matrix.
T make_shared(T... args)

◆ operator/()

std::shared_ptr< BoundedField > lsst::afw::math::BoundedField::operator/ ( double scale) const
inlineinherited

Definition at line 201 of file BoundedField.h.

201{ return (*this) * (1.0 / scale); }
scale(algorithm, min, max=None, frame=None)
Definition ds9.py:108

◆ operator=() [1/2]

ChebyshevBoundedField & lsst::afw::math::ChebyshevBoundedField::operator= ( ChebyshevBoundedField && )
delete

◆ operator=() [2/2]

ChebyshevBoundedField & lsst::afw::math::ChebyshevBoundedField::operator= ( ChebyshevBoundedField const & )
delete

◆ operator==()

bool lsst::afw::math::ChebyshevBoundedField::operator== ( BoundedField const & rhs) const
overridevirtual

BoundedFields (of the same sublcass) are equal if their bounding boxes and parameters are equal.

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

Definition at line 394 of file ChebyshevBoundedField.cc.

394 {
395 auto rhsCasted = dynamic_cast<ChebyshevBoundedField const*>(&rhs);
396 if (!rhsCasted) return false;
397
398 return (getBBox() == rhsCasted->getBBox()) &&
399 (_coefficients.getShape() == rhsCasted->_coefficients.getShape()) &&
400 all(equal(_coefficients, rhsCasted->_coefficients));
401}
T equal(T... args)
bool all(CoordinateExpr< N > const &expr) noexcept
Return true if all elements are true.

◆ readFits() [1/6]

static std::shared_ptr< BoundedField > lsst::afw::table::io::PersistableFacade< BoundedField >::readFits ( fits::Fits & fitsfile)
inlinestaticinherited

Read an object from an already open FITS object.

Parameters
[in]fitsfileFITS object to read from, already positioned at the desired HDU.

Definition at line 183 of file Persistable.h.

183 {
185 }
static std::shared_ptr< BoundedField > dynamicCast(std::shared_ptr< Persistable > const &ptr)

◆ readFits() [2/6]

static std::shared_ptr< BoundedField > lsst::afw::table::io::PersistableFacade< BoundedField >::readFits ( fits::MemFileManager & manager,
int hdu = fits::DEFAULT_HDU )
inlinestaticinherited

Read an object from a FITS file in memory.

Parameters
[in]managerManager for the memory to read from.
[in]hduHDU to read, where 0 is the primary. The special value of afw::fits::DEFAULT_HDU skips the primary HDU if it is empty.

Definition at line 205 of file Persistable.h.

205 {
207 }

◆ readFits() [3/6]

static std::shared_ptr< BoundedField > lsst::afw::table::io::PersistableFacade< BoundedField >::readFits ( std::string const & fileName,
int hdu = fits::DEFAULT_HDU )
inlinestaticinherited

Read an object from a regular FITS file.

Parameters
[in]fileNameName of the file to read.
[in]hduHDU to read, where 0 is the primary. The special value of afw::fits::DEFAULT_HDU skips the primary HDU if it is empty.

Definition at line 194 of file Persistable.h.

194 {
196 }

◆ readFits() [4/6]

Read an object from an already open FITS object.

Parameters
[in]fitsfileFITS object to read from, already positioned at the desired HDU.

Definition at line 183 of file Persistable.h.

183 {
185 }

◆ readFits() [5/6]

Read an object from a FITS file in memory.

Parameters
[in]managerManager for the memory to read from.
[in]hduHDU to read, where 0 is the primary. The special value of afw::fits::DEFAULT_HDU skips the primary HDU if it is empty.

Definition at line 205 of file Persistable.h.

205 {
207 }

◆ readFits() [6/6]

static std::shared_ptr< ChebyshevBoundedField > lsst::afw::table::io::PersistableFacade< ChebyshevBoundedField >::readFits ( std::string const & fileName,
int hdu = fits::DEFAULT_HDU )
inlinestaticinherited

Read an object from a regular FITS file.

Parameters
[in]fileNameName of the file to read.
[in]hduHDU to read, where 0 is the primary. The special value of afw::fits::DEFAULT_HDU skips the primary HDU if it is empty.

Definition at line 194 of file Persistable.h.

194 {
196 }

◆ relocate()

std::shared_ptr< ChebyshevBoundedField > lsst::afw::math::ChebyshevBoundedField::relocate ( lsst::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 254 of file ChebyshevBoundedField.cc.

254 {
255 return std::make_shared<ChebyshevBoundedField>(bbox, _coefficients);
256}

◆ toString()

std::string lsst::afw::math::ChebyshevBoundedField::toString ( ) const
overrideprivatevirtual

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

Definition at line 403 of file ChebyshevBoundedField.cc.

403 {
404 std::ostringstream os;
405 os << "ChebyshevBoundedField (" << _coefficients.getShape() << " coefficients in y,x)";
406 return os.str();
407}
T str(T... args)

◆ truncate()

std::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 230 of file ChebyshevBoundedField.cc.

230 {
231 if (static_cast<std::size_t>(ctrl.orderX) >= _coefficients.getSize<1>()) {
232 throw LSST_EXCEPT(pex::exceptions::LengthError,
233 (boost::format("New x order (%d) exceeds old x order (%d)") % ctrl.orderX %
234 (_coefficients.getSize<1>() - 1))
235 .str());
236 }
237 if (static_cast<std::size_t>(ctrl.orderY) >= _coefficients.getSize<0>()) {
238 throw LSST_EXCEPT(pex::exceptions::LengthError,
239 (boost::format("New y order (%d) exceeds old y order (%d)") % ctrl.orderY %
240 (_coefficients.getSize<0>() - 1))
241 .str());
242 }
243 ndarray::Array<double, 2, 2> coefficients = ndarray::allocate(ctrl.orderY + 1, ctrl.orderX + 1);
244 coefficients.deep() = _coefficients[ndarray::view(0, ctrl.orderY + 1)(0, ctrl.orderX + 1)];
245 if (ctrl.triangular) {
246 Packer packer(ctrl);
247 ndarray::Array<double, 1, 1> packed = ndarray::allocate(packer.size);
248 packer.pack(packed, coefficients);
249 packer.unpack(coefficients, packed);
250 }
251 return std::make_shared<ChebyshevBoundedField>(getBBox(), coefficients);
252}
decltype(sizeof(void *)) size_t
Definition doctest.h:524

◆ write()

void lsst::afw::math::ChebyshevBoundedField::write ( OutputArchiveHandle & handle) const
overrideprotectedvirtual

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 378 of file ChebyshevBoundedField.cc.

378 {
379 PersistenceHelper const keys(_coefficients.getSize<1>(), _coefficients.getSize<0>());
380 table::BaseCatalog catalog = handle.makeCatalog(keys.schema);
381 std::shared_ptr<table::BaseRecord> record = catalog.addNew();
382 record->set(keys.orderX, _coefficients.getSize<1>() - 1);
383 record->set(keys.bbox, getBBox());
384 (*record)[keys.coefficients].deep() = ndarray::flatten<1>(_coefficients);
385 handle.saveCatalog(catalog);
386}
CatalogT< BaseRecord > BaseCatalog
Definition fwd.h:72

◆ writeFits() [1/3]

void lsst::afw::table::io::Persistable::writeFits ( fits::Fits & fitsfile) const
inherited

Write the object to an already-open FITS object.

Parameters
[in]fitsfileOpen FITS object to write to.

Definition at line 18 of file Persistable.cc.

18 {
19 OutputArchive archive;
20 archive.put(this);
21 archive.writeFits(fitsfile);
22}

◆ writeFits() [2/3]

void lsst::afw::table::io::Persistable::writeFits ( fits::MemFileManager & manager,
std::string const & mode = "w" ) const
inherited

Write the object to a FITS image in memory.

Parameters
[in]managerName of the file to write to.
[in]modeIf "w", any existing file with the given name will be overwritten. If "a", new HDUs will be appended to an existing file.

Definition at line 29 of file Persistable.cc.

29 {
30 fits::Fits fitsfile(manager, mode, fits::Fits::AUTO_CLOSE | fits::Fits::AUTO_CHECK);
31 writeFits(fitsfile);
32}
void writeFits(std::string const &fileName, std::string const &mode="w") const
Write the object to a regular FITS file.

◆ writeFits() [3/3]

void lsst::afw::table::io::Persistable::writeFits ( std::string const & fileName,
std::string const & mode = "w" ) const
inherited

Write the object to a regular FITS file.

Parameters
[in]fileNameName of the file to write to.
[in]modeIf "w", any existing file with the given name will be overwritten. If "a", new HDUs will be appended to an existing file.

Definition at line 24 of file Persistable.cc.

24 {
25 fits::Fits fitsfile(fileName, mode, fits::Fits::AUTO_CLOSE | fits::Fits::AUTO_CHECK);
26 writeFits(fitsfile);
27}

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