LSST Applications g0b6bd0c080+a72a5dd7e6,g1182afd7b4+2a019aa3bb,g17e5ecfddb+2b8207f7de,g1d67935e3f+06cf436103,g38293774b4+ac198e9f13,g396055baef+6a2097e274,g3b44f30a73+6611e0205b,g480783c3b1+98f8679e14,g48ccf36440+89c08d0516,g4b93dc025c+98f8679e14,g5c4744a4d9+a302e8c7f0,g613e996a0d+e1c447f2e0,g6c8d09e9e7+25247a063c,g7271f0639c+98f8679e14,g7a9cd813b8+124095ede6,g9d27549199+a302e8c7f0,ga1cf026fa3+ac198e9f13,ga32aa97882+7403ac30ac,ga786bb30fb+7a139211af,gaa63f70f4e+9994eb9896,gabf319e997+ade567573c,gba47b54d5d+94dc90c3ea,gbec6a3398f+06cf436103,gc6308e37c7+07dd123edb,gc655b1545f+ade567573c,gcc9029db3c+ab229f5caf,gd01420fc67+06cf436103,gd877ba84e5+06cf436103,gdb4cecd868+6f279b5b48,ge2d134c3d5+cc4dbb2e3f,ge448b5faa6+86d1ceac1d,gecc7e12556+98f8679e14,gf3ee170dca+25247a063c,gf4ac96e456+ade567573c,gf9f5ea5b4d+ac198e9f13,gff490e6085+8c2580be5c,w.2022.27
LSST Data Management Base Package
ConstrainedPhotometryModel.h
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#ifndef LSST_JOINTCAL_CONSTRAINED_PHOTOMETRY_MODEL_H
26#define LSST_JOINTCAL_CONSTRAINED_PHOTOMETRY_MODEL_H
27
28#include <map>
29#include <utility>
30
34
35namespace lsst {
36namespace jointcal {
37
50public:
61 explicit ConstrainedPhotometryModel(CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox,
62 LOG_LOGGER log, int visitOrder = 7, double errorPedestal = 0)
63 : PhotometryModel(std::move(log), errorPedestal), _fittingChips(false), _fittingVisits(false) {
64 _chipVisitMap.reserve(ccdImageList.size());
65 }
66
72
74 Eigen::Index assignIndices(std::string const &whatToFit, Eigen::Index firstIndex) override;
75
77 void offsetParams(Eigen::VectorXd const &delta) override;
78
80 void freezeErrorTransform() override;
81
83 void getMappingIndices(CcdImage const &ccdImage, IndexVector &indices) const override;
84
86 std::size_t getTotalParameters() const override;
87
89 void computeParameterDerivatives(MeasuredStar const &measuredStar, CcdImage const &ccdImage,
90 Eigen::VectorXd &derivatives) const override;
91
93 void print(std::ostream &out) const override;
94
96
97protected:
98 PhotometryMappingBase *findMapping(CcdImage const &ccdImage) const override;
99
100 /* The per-ccdImage transforms, each of which is a composition of a chip and visit transform.
101 * Not all pairs of _visitMap[visit] and _chipMap[chip] are guaranteed to have an entry in
102 * _chipVisitMap (for example, one ccd in one visit might fail to produce a catalog).
103 */
106
107 // The per-visit transforms that go into _chipVisitMap.
110 // The per-chip transforms that go into _chipVisitMap.
113
117 template <class ChipTransform, class VisitTransform, class ChipVisitMapping>
118 void initialize(CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox, int visitOrder);
119
122
128 double visitMean;
129 };
130
132 PrepPhotoCalib prepPhotoCalib(CcdImage const &ccdImage) const;
133
134private:
135 // Which components of the model are we fitting currently?
136 bool _fittingChips;
137 bool _fittingVisits;
138};
139
141public:
142 explicit ConstrainedFluxModel(CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox,
143 int visitOrder = 7, double errorPedestal = 0)
144 : ConstrainedPhotometryModel(ccdImageList, focalPlaneBBox,
145 LOG_GET("lsst.jointcal.ConstrainedFluxModel"), visitOrder,
147 initialize<FluxTransformSpatiallyInvariant, FluxTransformChebyshev, ChipVisitFluxMapping>(
148 ccdImageList, focalPlaneBBox, visitOrder);
149 }
150
152 void offsetFittedStar(FittedStar &fittedStar, double delta) const override {
153 fittedStar.getFlux() -= delta;
154 }
155
157 double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
158
160 double computeRefResidual(FittedStar const &fittedStar, RefStar const &refStar) const override {
161 return fittedStar.getFlux() - refStar.getFlux();
162 };
163
165 double getRefError(RefStar const &refStar) const override { return refStar.getFluxErr(); }
166
168 double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
169
171 double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
172
174 std::shared_ptr<afw::image::PhotoCalib> toPhotoCalib(CcdImage const &ccdImage) const override;
175
177 void print(std::ostream &out) const override;
178
179protected:
182 return photoCalib->getCalibrationMean();
183 }
184};
185
187public:
188 explicit ConstrainedMagnitudeModel(CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox,
189 int visitOrder = 7, double errorPedestal = 0)
190 : ConstrainedPhotometryModel(ccdImageList, focalPlaneBBox,
191 LOG_GET("lsst.jointcal.ConstrainedMagnitudeModel"), visitOrder,
194 ChipVisitMagnitudeMapping>(ccdImageList, focalPlaneBBox, visitOrder);
195 }
196
198 void offsetFittedStar(FittedStar &fittedStar, double delta) const override {
199 fittedStar.getMag() -= delta;
200 }
201
203 double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
204
206 double computeRefResidual(FittedStar const &fittedStar, RefStar const &refStar) const override {
207 return fittedStar.getMag() - refStar.getMag();
208 };
209
211 double getRefError(RefStar const &refStar) const override { return refStar.getMagErr(); }
212
214 double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
215
217 double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
218
220 std::shared_ptr<afw::image::PhotoCalib> toPhotoCalib(CcdImage const &ccdImage) const override;
221
223 void print(std::ostream &out) const override;
224
225protected:
228 return utils::nanojanskyToABMagnitude(photoCalib->getCalibrationMean());
229 }
230};
231
232} // namespace jointcal
233} // namespace lsst
234
235#endif // LSST_JOINTCAL_CONSTRAINED_PHOTOMETRY_MODEL_H
#define LOG_GET(logger)
Returns a Log object associated with logger.
Definition: Log.h:75
#define LOG_LOGGER
Definition: Log.h:714
Transform LSST spatial data, such as lsst::geom::Point2D and lsst::geom::SpherePoint,...
Definition: Transform.h:68
A floating-point coordinate rectangle geometry.
Definition: Box.h:413
double getMagErr() const
Definition: BaseStar.h:108
double getMag() const
Definition: BaseStar.h:105
double getFlux() const
Definition: BaseStar.h:98
double getFluxErr() const
Definition: BaseStar.h:102
Handler of an actual image from a single CCD.
Definition: CcdImage.h:64
std::shared_ptr< afw::image::PhotoCalib > toPhotoCalib(CcdImage const &ccdImage) const override
Return the mapping of ccdImage represented as a PhotoCalib.
double initialChipCalibration(std::shared_ptr< afw::image::PhotoCalib const > photoCalib) override
Return the initial calibration to use from this photoCalib.
double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Compute the residual between the model applied to a star and its associated fittedStar.
ConstrainedFluxModel(CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox, int visitOrder=7, double errorPedestal=0)
void offsetFittedStar(FittedStar &fittedStar, double delta) const override
Offset the appropriate flux or magnitude (by -delta).
double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux uncertainty for measuredStar on ccdImage.
double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux for measuredStar on ccdImage.
void print(std::ostream &out) const override
Print a string representation of the contents of this mapping, for debugging.
double getRefError(RefStar const &refStar) const override
Return the refStar error appropriate for this model (e.g. fluxErr or magErr).
double computeRefResidual(FittedStar const &fittedStar, RefStar const &refStar) const override
Return the fittedStar - refStar residual appropriate for this model (e.g. flux - flux or mag - mag).
double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Compute the residual between the model applied to a star and its associated fittedStar.
void offsetFittedStar(FittedStar &fittedStar, double delta) const override
Offset the appropriate flux or magnitude (by -delta).
ConstrainedMagnitudeModel(CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox, int visitOrder=7, double errorPedestal=0)
double computeRefResidual(FittedStar const &fittedStar, RefStar const &refStar) const override
Return the fittedStar - refStar residual appropriate for this model (e.g. flux - flux or mag - mag).
double initialChipCalibration(std::shared_ptr< afw::image::PhotoCalib const > photoCalib) override
Return the initial calibration to use from this photoCalib.
double getRefError(RefStar const &refStar) const override
Return the refStar error appropriate for this model (e.g. fluxErr or magErr).
void print(std::ostream &out) const override
Print a string representation of the contents of this mapping, for debugging.
double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux uncertainty for measuredStar on ccdImage.
double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux for measuredStar on ccdImage.
std::shared_ptr< afw::image::PhotoCalib > toPhotoCalib(CcdImage const &ccdImage) const override
Return the mapping of ccdImage represented as a PhotoCalib.
PrepPhotoCalib prepPhotoCalib(CcdImage const &ccdImage) const
Helper for preparing toPhotoCalib()
void initialize(CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox, int visitOrder)
Initialize the chip, visit, and chipVisit mappings by creating appropriate transforms and mappings.
PhotometryMappingBase * findMapping(CcdImage const &ccdImage) const override
Return a pointer to the mapping associated with this ccdImage.
void freezeErrorTransform() override
Once this routine has been called, the error transform is not modified by offsetParams().
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the provided amounts (by -delta).
ConstrainedPhotometryModel(ConstrainedPhotometryModel &&)=delete
ConstrainedPhotometryModel & operator=(ConstrainedPhotometryModel &&)=delete
ConstrainedPhotometryModel(CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox, LOG_LOGGER log, int visitOrder=7, double errorPedestal=0)
Construct a constrained photometry model.
void print(std::ostream &out) const override
Print a string representation of the contents of this mapping, for debugging.
std::size_t getTotalParameters() const override
Return the total number of parameters in this model.
Eigen::Index assignIndices(std::string const &whatToFit, Eigen::Index firstIndex) override
Assign indices in the full matrix to the parameters being fit in the mappings, starting at firstIndex...
void getMappingIndices(CcdImage const &ccdImage, IndexVector &indices) const override
Get how this set of parameters (of length Npar()) map into the "grand" fit.
ConstrainedPhotometryModel & operator=(ConstrainedPhotometryModel const &)=delete
void computeParameterDerivatives(MeasuredStar const &measuredStar, CcdImage const &ccdImage, Eigen::VectorXd &derivatives) const override
Compute the parametric derivatives of this model.
ConstrainedPhotometryModel(ConstrainedPhotometryModel const &)=delete
No copy or move: there is only ever one instance of a given model (i.e. per ccd+visit)
virtual double initialChipCalibration(std::shared_ptr< afw::image::PhotoCalib const > photoCalib)=0
Return the initial calibration to use from this photoCalib.
FittedStars are objects whose position or flux is going to be fitted, and which come from the associa...
Definition: FittedStar.h:50
nth-order 2d Chebyshev photometry transform, plus the input flux.
Photometric offset independent of position, defined as -2.5 * log(flux / fluxMag0).
Sources measured on images.
Definition: MeasuredStar.h:51
Relates transform(s) to their position in the fitting matrix and allows interaction with the transfor...
Objects used as position/flux anchors (e.g.
Definition: RefStar.h:45
T move(T... args)
double nanojanskyToABMagnitude(double flux)
Convert a flux in nanojansky to AB magnitude.
Definition: Magnitude.cc:30
A base class for image defects.
STL namespace.
T reserve(T... args)
T size(T... args)
std::shared_ptr< afw::geom::TransformPoint2ToPoint2 > pixToFocal
Key< int > photoCalib
Definition: Exposure.cc:67