LSST Applications g180d380827+770a9040cc,g2079a07aa2+86d27d4dc4,g2305ad1205+09cfdadad9,g2bbee38e9b+c6a8a0fb72,g337abbeb29+c6a8a0fb72,g33d1c0ed96+c6a8a0fb72,g3a166c0a6a+c6a8a0fb72,g3ddfee87b4+1ea5e09c42,g48712c4677+7e2ea9cd42,g487adcacf7+301d09421d,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+96fcb956a6,g64a986408d+23540ee355,g858d7b2824+23540ee355,g864b0138d7+aa38e45daa,g95921f966b+d83dc58ecd,g991b906543+23540ee355,g99cad8db69+7f13b58a93,g9c22b2923f+e2510deafe,g9ddcbc5298+9a081db1e4,ga1e77700b3+03d07e1c1f,gb0e22166c9+60f28cb32d,gb23b769143+23540ee355,gba4ed39666+c2a2e4ac27,gbb8dafda3b+49e7449578,gbd998247f1+585e252eca,gc120e1dc64+1bbfa184e1,gc28159a63d+c6a8a0fb72,gc3e9b769f7+385ea95214,gcf0d15dbbd+1ea5e09c42,gdaeeff99f8+f9a426f77a,ge6526c86ff+1bccc98490,ge79ae78c31+c6a8a0fb72,gee10cc3b42+585e252eca,w.2024.18
LSST Data Management Base Package
Loading...
Searching...
No Matches
PhotometryTransform.cc
Go to the documentation of this file.
1// -*- LSST-C++ -*-
2/*
3 * This file is part of jointcal.
4 *
5 * Developed for the LSST Data Management System.
6 * This product includes software developed by the LSST Project
7 * (https://www.lsst.org).
8 * See the COPYRIGHT file at the top-level directory of this distribution
9 * for details of code ownership.
10 *
11 * This program is free software: you can redistribute it and/or modify
12 * it under the terms of the GNU General Public License as published by
13 * the Free Software Foundation, either version 3 of the License, or
14 * (at your option) any later version.
15 *
16 * This program is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 * GNU General Public License for more details.
20 *
21 * You should have received a copy of the GNU General Public License
22 * along with this program. If not, see <https://www.gnu.org/licenses/>.
23 */
24
25#include "ndarray.h"
26#include "Eigen/Core"
27
29
30#include "lsst/jointcal/Point.h"
32
33namespace lsst {
34namespace jointcal {
35
36// ------------------ PhotometryTransformChebyshev helpers ---------------------------------------------------
37
38namespace {
39
40// To evaluate a 1-d Chebyshev function without needing to have workspace, we use the
41// Clenshaw algorith, which is like going through the recurrence relation in reverse.
42// The CoeffGetter argument g is something that behaves like an array, providing access
43// to the coefficients.
44template <typename CoeffGetter>
45double evaluateFunction1d(CoeffGetter g, double x, std::size_t size) {
46 double b_kp2 = 0.0, b_kp1 = 0.0;
47 for (std::size_t k = (size - 1); k > 0; --k) {
48 double b_k = g[k] + 2 * x * b_kp1 - b_kp2;
49 b_kp2 = b_kp1;
50 b_kp1 = b_k;
51 }
52 return g[0] + x * b_kp1 - b_kp2;
53}
54
55// This class imitates a 1-d array, by running evaluateFunction1d on a nested dimension;
56// this lets us reuse the logic in evaluateFunction1d for both dimensions. Essentially,
57// we run evaluateFunction1d on a column of coefficients to evaluate T_i(x), then pass
58// the result of that to evaluateFunction1d with the results as the "coefficients" associated
59// with the T_j(y) functions.
60struct RecursionArrayImitator {
61 double operator[](Eigen::Index i) const {
62 return evaluateFunction1d(coefficients[i], x, coefficients.getSize<1>());
63 }
64
65 RecursionArrayImitator(ndarray::Array<double const, 2, 2> const &coefficients_, double x_)
66 : coefficients(coefficients_), x(x_) {}
67
68 ndarray::Array<double const, 2, 2> coefficients;
69 double x;
70};
71
72// Compute an affine transform that maps an arbitrary box to [-1,1]x[-1,1]
73geom::AffineTransform makeChebyshevRangeTransform(geom::Box2D const &bbox) {
75 geom::LinearTransform::makeScaling(2.0 / bbox.getWidth(), 2.0 / bbox.getHeight()),
76 geom::Extent2D(-(2.0 * bbox.getCenterX()) / bbox.getWidth(),
77 -(2.0 * bbox.getCenterY()) / bbox.getHeight()));
78}
79
80// Initialize a "unit" Chebyshev
81ndarray::Array<double, 2, 2> _initializeChebyshev(size_t order, bool identity) {
82 ndarray::Array<double, 2, 2> coeffs = ndarray::allocate(ndarray::makeVector(order + 1, order + 1));
83 coeffs.deep() = 0.0;
84 if (identity) {
85 coeffs[0][0] = 1;
86 }
87 return coeffs;
88}
89} // namespace
90
92 bool identity)
93 : _bbox(bbox),
94 _toChebyshevRange(makeChebyshevRangeTransform(bbox)),
95 _coefficients(_initializeChebyshev(order, identity)),
96 _order(order),
97 _nParameters((order + 1) * (order + 2) / 2) {}
98
100 geom::Box2D const &bbox)
101 : _bbox(bbox),
102 _toChebyshevRange(makeChebyshevRangeTransform(bbox)),
103 _coefficients(coefficients),
104 _order(coefficients.size() - 1),
105 _nParameters((_order + 1) * (_order + 2) / 2) {}
106
107void PhotometryTransformChebyshev::offsetParams(Eigen::VectorXd const &delta) {
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}
117
118namespace {
119// The integral of T_n(x) over [-1,1]:
120// https://en.wikipedia.org/wiki/Chebyshev_polynomials#Differentiation_and_integration
121double integrateTn(int n) {
122 if (n % 2 == 1)
123 return 0;
124 else
125 return 2.0 / (1.0 - static_cast<double>(n * n));
126}
127
136Eigen::VectorXd computeIntegralTn(ndarray::Size order, double x) {
137 // We need to compute through order+2, because the integral recurrance relation uses Tn+1(x).
138 Eigen::VectorXd Tn(order + 2);
139
140 // The nth chebyshev polynomial evaluated at x.
141 Tn[0] = 1;
142 Tn[1] = x;
143 for (ndarray::Size i = 2; i <= order + 1; ++i) {
144 Tn[i] = 2 * x * Tn[i - 1] - Tn[i - 2];
145 }
146
147 // The integral of the nth chebyshev polynomial at x.
148 Eigen::VectorXd integralTn(order + 1);
149 integralTn[0] = x;
150 if (order > 0) { // have to evaluate this separately from the recurrence relation below.
151 integralTn[1] = 0.5 * x * x;
152 }
153 for (ndarray::Size i = 2; i <= order; ++i) {
154 // recurrence relation: integral(Tn(x)) = n*T[n+1](x)/(n^2 - 1) - x*T[n](x)/(n-1)
155 integralTn[i] = i * Tn[i + 1] / (i * i - 1) - x * Tn[i] / (i - 1);
156 }
157
158 return integralTn;
159}
160
161} // namespace
162
163double PhotometryTransformChebyshev::oneIntegral(double x, double y) const {
164 geom::Point2D point = _toChebyshevRange(geom::Point2D(x, y));
165
166 auto integralTnx = computeIntegralTn(_order, point.getX());
167 auto integralTmy = computeIntegralTn(_order, point.getY());
168
169 // NOTE: the indexing in this method and offsetParams must be kept consistent!
170 // Roll up the x and y terms
171 double result = 0;
172 for (ndarray::Size j = 0; j <= _order; ++j) {
173 ndarray::Size const iMax = _order - j; // to save re-computing `i+j <= order` every inner step.
174 for (ndarray::Size i = 0; i <= iMax; ++i) {
175 result += _coefficients[j][i] * integralTnx[i] * integralTmy[j];
176 }
177 }
178 return result;
179}
180
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}
192
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}
203
205 return integrate(bbox) / bbox.getArea();
206}
207
208double PhotometryTransformChebyshev::mean() const { return integrate() / _bbox.getArea(); }
209
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}
223
225 geom::Point2D p = _toChebyshevRange(geom::Point2D(x, y));
226 return evaluateFunction1d(RecursionArrayImitator(_coefficients, p.getX()), p.getY(),
227 _coefficients.getSize<0>());
228}
229
231 double x, double y, Eigen::Ref<Eigen::VectorXd> derivatives) const {
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}
257
258} // namespace jointcal
259} // namespace lsst
py::object result
Definition _schema.cc:429
AmpInfoBoxKey bbox
Definition Amplifier.cc:117
ndarray::Array< double const, 2, 2 > coefficients
int y
Definition SpanSet.cc:48
An affine coordinate transformation consisting of a linear transformation and an offset.
LinearTransform const & getLinear() const noexcept
A floating-point coordinate rectangle geometry.
Definition Box.h:413
double getArea() const noexcept
1-d interval accessors
Definition Box.h:531
static LinearTransform makeScaling(double s) noexcept
double computeDeterminant() const noexcept
Return the determinant of the 2x2 matrix.
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by some (negative) amount during fitting.
Eigen::VectorXd getParameters() const override
Get a copy of the parameters of this model, in the same order as offsetParams.
PhotometryTransformChebyshev(size_t order, geom::Box2D const &bbox, bool identity)
Create a Chebyshev transform with terms up to order in (x*y).
void computeChebyshevDerivatives(double x, double y, Eigen::Ref< Eigen::VectorXd > derivatives) const
Set the derivatives of this polynomial at x,y.
double computeChebyshev(double x, double y) const
Return the value of this polynomial at x,y.
table::Key< int > order