LSST Applications g0265f82a02+0e5473021a,g02d81e74bb+f5613e8b4f,g1470d8bcf6+190ad2ba91,g14a832a312+311607e4ab,g2079a07aa2+86d27d4dc4,g2305ad1205+a8e3196225,g295015adf3+b67ee847e5,g2bbee38e9b+0e5473021a,g337abbeb29+0e5473021a,g3ddfee87b4+a761f810f3,g487adcacf7+17c8fdbcbd,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+65b5bd823e,g5a732f18d5+53520f316c,g64a986408d+f5613e8b4f,g6c1bc301e9+51106c2951,g858d7b2824+f5613e8b4f,g8a8a8dda67+585e252eca,g99cad8db69+6729933424,g9ddcbc5298+9a081db1e4,ga1e77700b3+15fc3df1f7,ga8c6da7877+ef4e3a5875,gb0e22166c9+60f28cb32d,gb6a65358fc+0e5473021a,gba4ed39666+c2a2e4ac27,gbb8dafda3b+e9bba80f27,gc120e1dc64+eee469a5e5,gc28159a63d+0e5473021a,gcf0d15dbbd+a761f810f3,gdaeeff99f8+f9a426f77a,ge6526c86ff+d4c1d4bfef,ge79ae78c31+0e5473021a,gee10cc3b42+585e252eca,gf1cff7945b+f5613e8b4f,w.2024.16
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | List of all members
lsst::jointcal::PhotometryTransformChebyshev Class Referenceabstract

nth-order 2d Chebyshev photometry transform. More...

#include <PhotometryTransform.h>

Inheritance diagram for lsst::jointcal::PhotometryTransformChebyshev:
lsst::jointcal::PhotometryTransform lsst::jointcal::FluxTransformChebyshev lsst::jointcal::MagnitudeTransformChebyshev

Public Member Functions

 PhotometryTransformChebyshev (size_t order, geom::Box2D const &bbox, bool identity)
 Create a Chebyshev transform with terms up to order in (x*y).
 
 PhotometryTransformChebyshev (ndarray::Array< double, 2, 2 > const &coefficients, geom::Box2D const &bbox)
 Create a Chebyshev transform with the specified coefficients.
 
double transformError (double x, double y, double value, double valueErr) const override
 Return the transformed valueErr at Point(x,y).
 
void print (std::ostream &out) const override
 Print the transform coefficients to stream.
 
std::size_t getNpar () const override
 Return the number of parameters (used to compute chisq)
 
void offsetParams (Eigen::VectorXd const &delta) override
 Offset the parameters by some (negative) amount during fitting.
 
ndarray::Array< double, 2, 2 > getCoefficients () const
 Get a copy of the coefficients of the polynomials, as a 2d array (NOTE: layout is [y][x])
 
Eigen::VectorXd getParameters () const override
 Get a copy of the parameters of this model, in the same order as offsetParams.
 
ndarray::Size getOrder () const
 
geom::Box2D getBBox () const
 
double mean (geom::Box2D const &bbox) const
 Compute the mean of this tranform on the bbox (default to our bbox).
 
double mean () const
 
double integrate (geom::Box2D const &bbox) const
 
double integrate () const
 
virtual double transform (double x, double y, double value) const =0
 Return the transform of value at (x,y).
 
double transform (Point const &in, double value) const
 Return the transformed value at Point(x,y).
 
double transformError (Point const &in, double value, double valueErr) const
 Return the transformed valueErr at Point(x,y).
 
virtual std::shared_ptr< PhotometryTransformclone () const =0
 return a copy (allocated by new) of the transformation.
 
virtual void computeParameterDerivatives (double x, double y, double value, Eigen::Ref< Eigen::VectorXd > derivatives) const =0
 Compute the derivatives with respect to the parameters (i.e.
 

Protected Member Functions

double computeChebyshev (double x, double y) const
 Return the value of this polynomial at x,y.
 
void computeChebyshevDerivatives (double x, double y, Eigen::Ref< Eigen::VectorXd > derivatives) const
 Set the derivatives of this polynomial at x,y.
 

Detailed Description

nth-order 2d Chebyshev photometry transform.

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.

Note that the polynomial order=n means that the highest terms will be of the form:

\[ a_{0,n}*x^n*y^0, a_{n-1,1}*x^(n-1)*y^1, ..., a_{1,n-1}*x^1*y^(n-1), a_{n,0}*x^0*y^n \]

Definition at line 222 of file PhotometryTransform.h.

Constructor & Destructor Documentation

◆ PhotometryTransformChebyshev() [1/2]

lsst::jointcal::PhotometryTransformChebyshev::PhotometryTransformChebyshev ( size_t order,
geom::Box2D const & bbox,
bool identity )

Create a Chebyshev transform with terms up to order in (x*y).

Parameters
[in]orderThe maximum order in (x*y).
[in]bboxThe bounding box it is valid within, to rescale it to [-1,1].
[in]identityIf true, set a_0,0==1, otherwise all coefficients are 0.

Definition at line 91 of file PhotometryTransform.cc.

93 : _bbox(bbox),
94 _toChebyshevRange(makeChebyshevRangeTransform(bbox)),
95 _coefficients(_initializeChebyshev(order, identity)),
96 _order(order),
97 _nParameters((order + 1) * (order + 2) / 2) {}
AmpInfoBoxKey bbox
Definition Amplifier.cc:117
table::Key< int > order

◆ PhotometryTransformChebyshev() [2/2]

lsst::jointcal::PhotometryTransformChebyshev::PhotometryTransformChebyshev ( ndarray::Array< double, 2, 2 > const & coefficients,
geom::Box2D const & bbox )

Create a Chebyshev transform with the specified coefficients.

The polynomial order is determined from the number of coefficients, taking only the anti-diagonal upper triangle portion of the passed-in coefficients

Parameters
coefficientsThe polynomial coefficients.
[in]bboxThe bounding box it is valid within, to rescale it to [-1,1].

Definition at line 99 of file PhotometryTransform.cc.

101 : _bbox(bbox),
102 _toChebyshevRange(makeChebyshevRangeTransform(bbox)),
103 _coefficients(coefficients),
104 _order(coefficients.size() - 1),
105 _nParameters((_order + 1) * (_order + 2) / 2) {}
ndarray::Array< double const, 2, 2 > coefficients

Member Function Documentation

◆ clone()

virtual std::shared_ptr< PhotometryTransform > lsst::jointcal::PhotometryTransform::clone ( ) const
pure virtualinherited

◆ computeChebyshev()

double lsst::jointcal::PhotometryTransformChebyshev::computeChebyshev ( double x,
double y ) const
protected

Return the value of this polynomial at x,y.

For use in the sublcass transform() methods.

Definition at line 224 of file PhotometryTransform.cc.

224 {
225 geom::Point2D p = _toChebyshevRange(geom::Point2D(x, y));
226 return evaluateFunction1d(RecursionArrayImitator(_coefficients, p.getX()), p.getY(),
227 _coefficients.getSize<0>());
228}
int y
Definition SpanSet.cc:48

◆ computeChebyshevDerivatives()

void lsst::jointcal::PhotometryTransformChebyshev::computeChebyshevDerivatives ( double x,
double y,
Eigen::Ref< Eigen::VectorXd > derivatives ) const
protected

Set the derivatives of this polynomial at x,y.

For use in the sublcass computeParameterDerivatives() methods.

Definition at line 230 of file PhotometryTransform.cc.

231 {
232 geom::Point2D p = _toChebyshevRange(geom::Point2D(x, y));
233 // Algorithm: compute all the individual components recursively (since we'll need them anyway),
234 // then combine them into the final answer vectors.
235 Eigen::VectorXd Tnx(_order + 1);
236 Eigen::VectorXd Tmy(_order + 1);
237 Tnx[0] = 1;
238 Tmy[0] = 1;
239 if (_order >= 1) {
240 Tnx[1] = p.getX();
241 Tmy[1] = p.getY();
242 }
243 for (ndarray::Size i = 2; i <= _order; ++i) {
244 Tnx[i] = 2 * p.getX() * Tnx[i - 1] - Tnx[i - 2];
245 Tmy[i] = 2 * p.getY() * Tmy[i - 1] - Tmy[i - 2];
246 }
247
248 // NOTE: the indexing in this method and offsetParams must be kept consistent!
249 Eigen::VectorXd::Index k = 0;
250 for (ndarray::Size j = 0; j <= _order; ++j) {
251 ndarray::Size const iMax = _order - j; // to save re-computing `i+j <= order` every inner step.
252 for (ndarray::Size i = 0; i <= iMax; ++i, ++k) {
253 derivatives[k] = Tmy[j] * Tnx[i];
254 }
255 }
256}

◆ computeParameterDerivatives()

virtual void lsst::jointcal::PhotometryTransform::computeParameterDerivatives ( double x,
double y,
double value,
Eigen::Ref< Eigen::VectorXd > derivatives ) const
pure virtualinherited

Compute the derivatives with respect to the parameters (i.e.

the coefficients).

Parameters
[in]xThe x coordinate to compute at (in the appropriate units for this transform).
[in]yThe y coordinate to compute at (in the appropriate units for this transform).
[in]valueThe instrument flux or magnitude to compute the derivative at.
[out]derivativesThe computed derivatives, in the same order as the deltas in offsetParams.

Implemented in lsst::jointcal::FluxTransformSpatiallyInvariant, lsst::jointcal::MagnitudeTransformSpatiallyInvariant, lsst::jointcal::FluxTransformChebyshev, and lsst::jointcal::MagnitudeTransformChebyshev.

◆ getBBox()

geom::Box2D lsst::jointcal::PhotometryTransformChebyshev::getBBox ( ) const
inline

Definition at line 264 of file PhotometryTransform.h.

264{ return _bbox; }

◆ getCoefficients()

ndarray::Array< double, 2, 2 > lsst::jointcal::PhotometryTransformChebyshev::getCoefficients ( ) const
inline

Get a copy of the coefficients of the polynomials, as a 2d array (NOTE: layout is [y][x])

Definition at line 257 of file PhotometryTransform.h.

257{ return ndarray::copy(_coefficients); }

◆ getNpar()

std::size_t lsst::jointcal::PhotometryTransformChebyshev::getNpar ( ) const
inlineoverridevirtual

Return the number of parameters (used to compute chisq)

Implements lsst::jointcal::PhotometryTransform.

Definition at line 251 of file PhotometryTransform.h.

251{ return _nParameters; }

◆ getOrder()

ndarray::Size lsst::jointcal::PhotometryTransformChebyshev::getOrder ( ) const
inline

Definition at line 262 of file PhotometryTransform.h.

262{ return _order; }

◆ getParameters()

Eigen::VectorXd lsst::jointcal::PhotometryTransformChebyshev::getParameters ( ) const
overridevirtual

Get a copy of the parameters of this model, in the same order as offsetParams.

Implements lsst::jointcal::PhotometryTransform.

Definition at line 210 of file PhotometryTransform.cc.

210 {
211 Eigen::VectorXd parameters(_nParameters);
212 // NOTE: the indexing in this method and offsetParams must be kept consistent!
213 Eigen::VectorXd::Index k = 0;
214 for (ndarray::Size j = 0; j <= _order; ++j) {
215 ndarray::Size const iMax = _order - j; // to save re-computing `i+j <= order` every inner step.
216 for (ndarray::Size i = 0; i <= iMax; ++i, ++k) {
217 parameters[k] = _coefficients[j][i];
218 }
219 }
220
221 return parameters;
222}

◆ integrate() [1/2]

double lsst::jointcal::PhotometryTransformChebyshev::integrate ( ) const

Definition at line 193 of file PhotometryTransform.cc.

193 {
194 double result = 0;
195 double determinant = _bbox.getArea() / 4.0;
196 for (ndarray::Size j = 0; j < _coefficients.getSize<0>(); j++) {
197 for (ndarray::Size i = 0; i < _coefficients.getSize<1>(); i++) {
198 result += _coefficients[j][i] * integrateTn(i) * integrateTn(j);
199 }
200 }
201 return result * determinant;
202}
py::object result
Definition _schema.cc:429
double getArea() const noexcept
1-d interval accessors
Definition Box.h:531

◆ integrate() [2/2]

double lsst::jointcal::PhotometryTransformChebyshev::integrate ( geom::Box2D const & bbox) const

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 181 of file PhotometryTransform.cc.

181 {
182 double result = 0;
183
184 result += oneIntegral(bbox.getMaxX(), bbox.getMaxY());
185 result += oneIntegral(bbox.getMinX(), bbox.getMinY());
186 result -= oneIntegral(bbox.getMaxX(), bbox.getMinY());
187 result -= oneIntegral(bbox.getMinX(), bbox.getMaxY());
188
189 // scale factor due to the change of limits in the integral
190 return result / _toChebyshevRange.getLinear().computeDeterminant();
191}
LinearTransform const & getLinear() const noexcept
double computeDeterminant() const noexcept
Return the determinant of the 2x2 matrix.

◆ mean() [1/2]

double lsst::jointcal::PhotometryTransformChebyshev::mean ( ) const

Definition at line 208 of file PhotometryTransform.cc.

208{ return integrate() / _bbox.getArea(); }

◆ mean() [2/2]

double lsst::jointcal::PhotometryTransformChebyshev::mean ( geom::Box2D const & bbox) const

Compute the mean of this tranform on the bbox (default to our bbox).

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 204 of file PhotometryTransform.cc.

204 {
205 return integrate(bbox) / bbox.getArea();
206}

◆ offsetParams()

void lsst::jointcal::PhotometryTransformChebyshev::offsetParams ( Eigen::VectorXd const & delta)
overridevirtual

Offset the parameters by some (negative) amount during fitting.

Equivalent to flatten(parameters) -= delta

Ordering of delta is the same as the ordering of the derivatives returned from computeParameterDerivatives.

Implements lsst::jointcal::PhotometryTransform.

Definition at line 107 of file PhotometryTransform.cc.

107 {
108 // NOTE: the indexing in this method and computeParameterDerivatives must be kept consistent!
109 Eigen::VectorXd::Index k = 0;
110 for (ndarray::Size j = 0; j <= _order; ++j) {
111 ndarray::Size const iMax = _order - j; // to save re-computing `i+j <= order` every inner step.
112 for (ndarray::Size i = 0; i <= iMax; ++i, ++k) {
113 _coefficients[j][i] -= delta[k];
114 }
115 }
116}

◆ print()

void lsst::jointcal::PhotometryTransformChebyshev::print ( std::ostream & out) const
inlineoverridevirtual

Print the transform coefficients to stream.

Implements lsst::jointcal::PhotometryTransform.

Definition at line 248 of file PhotometryTransform.h.

248{ out << "PhotometryTransformChebyshev: " << _coefficients; }

◆ transform() [1/2]

virtual double lsst::jointcal::PhotometryTransform::transform ( double x,
double y,
double value ) const
pure virtualinherited

◆ transform() [2/2]

double lsst::jointcal::PhotometryTransform::transform ( Point const & in,
double value ) const
inlineinherited

Return the transformed value at Point(x,y).

Definition at line 58 of file PhotometryTransform.h.

58{ return transform(in.x, in.y, value); }
table::Key< int > transform

◆ transformError() [1/2]

double lsst::jointcal::PhotometryTransformChebyshev::transformError ( double x,
double y,
double value,
double valueErr ) const
inlineoverridevirtual

Return the transformed valueErr at Point(x,y).

Implements lsst::jointcal::PhotometryTransform.

Definition at line 245 of file PhotometryTransform.h.

245{ return 0; }

◆ transformError() [2/2]

double lsst::jointcal::PhotometryTransform::transformError ( Point const & in,
double value,
double valueErr ) const
inlineinherited

Return the transformed valueErr at Point(x,y).

Definition at line 64 of file PhotometryTransform.h.

64 {
65 return transformError(in.x, in.y, value, valueErr);
66 }
virtual double transformError(double x, double y, double value, double valueErr) const =0
Return the transformed valueErr at Point(x,y).

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