LSSTApplications  18.0.0+106,18.0.0+50,19.0.0,19.0.0+1,19.0.0+10,19.0.0+11,19.0.0+13,19.0.0+17,19.0.0+2,19.0.0-1-g20d9b18+6,19.0.0-1-g425ff20,19.0.0-1-g5549ca4,19.0.0-1-g580fafe+6,19.0.0-1-g6fe20d0+1,19.0.0-1-g7011481+9,19.0.0-1-g8c57eb9+6,19.0.0-1-gb5175dc+11,19.0.0-1-gdc0e4a7+9,19.0.0-1-ge272bc4+6,19.0.0-1-ge3aa853,19.0.0-10-g448f008b,19.0.0-12-g6990b2c,19.0.0-2-g0d9f9cd+11,19.0.0-2-g3d9e4fb2+11,19.0.0-2-g5037de4,19.0.0-2-gb96a1c4+3,19.0.0-2-gd955cfd+15,19.0.0-3-g2d13df8,19.0.0-3-g6f3c7dc,19.0.0-4-g725f80e+11,19.0.0-4-ga671dab3b+1,19.0.0-4-gad373c5+3,19.0.0-5-ga2acb9c+2,19.0.0-5-gfe96e6c+2,w.2020.01
LSSTDataManagementBasePackage
SimplePhotometryModel.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_SIMPLE_PHOTOMETRY_MODEL_H
26 #define LSST_JOINTCAL_SIMPLE_PHOTOMETRY_MODEL_H
27 
28 #include <map>
29 
30 #include "lsst/jointcal/CcdImage.h"
34 #include "lsst/jointcal/Point.h"
35 
36 namespace lsst {
37 namespace jointcal {
38 
39 class CcdImage;
40 class Point;
41 
44 public:
45  SimplePhotometryModel(CcdImageList const &ccdImageList, LOG_LOGGER log, double errorPedestal_ = 0)
46  : PhotometryModel(log, errorPedestal_) {
47  _myMap.reserve(ccdImageList.size());
48  }
49 
55 
57  Eigen::Index assignIndices(std::string const &whatToFit, Eigen::Index firstIndex) override;
58 
60  void offsetParams(Eigen::VectorXd const &delta) override;
61 
63  void freezeErrorTransform() override;
64 
66  void getMappingIndices(CcdImage const &ccdImage, IndexVector &indices) const override;
67 
69  std::size_t getTotalParameters() const override;
70 
72  void computeParameterDerivatives(MeasuredStar const &measuredStar, CcdImage const &ccdImage,
73  Eigen::VectorXd &derivatives) const override;
74 
76  void dump(std::ostream &stream = std::cout) const override;
77 
78 protected:
80  MapType _myMap;
81 
83  PhotometryMappingBase *findMapping(CcdImage const &ccdImage) const override;
84 };
85 
87 public:
88  SimpleFluxModel(CcdImageList const &ccdImageList, double errorPedestal_ = 0);
89 
91  void offsetFittedStar(FittedStar &fittedStar, double delta) const override {
92  fittedStar.getFlux() -= delta;
93  }
94 
96  double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
97 
99  double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
100 
102  double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
103 
105  double getRefError(RefStar const &refStar) const override { return refStar.getFluxErr(); }
106 
108  double computeRefResidual(FittedStar const &fittedStar, RefStar const &refStar) const override {
109  return fittedStar.getFlux() - refStar.getFlux();
110  };
111 
117  std::shared_ptr<afw::image::PhotoCalib> toPhotoCalib(CcdImage const &ccdImage) const override;
118 };
119 
121 public:
122  SimpleMagnitudeModel(CcdImageList const &ccdImageList, double errorPedestal_ = 0);
123 
125  void offsetFittedStar(FittedStar &fittedStar, double delta) const override {
126  fittedStar.getMag() -= delta;
127  }
128 
130  double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
131 
133  double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
134 
136  double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override;
137 
139  double getRefError(RefStar const &refStar) const override { return refStar.getMagErr(); }
140 
142  double computeRefResidual(FittedStar const &fittedStar, RefStar const &refStar) const override {
143  return fittedStar.getMag() - refStar.getMag();
144  };
145 
151  std::shared_ptr<afw::image::PhotoCalib> toPhotoCalib(CcdImage const &ccdImage) const override;
152 };
153 
154 } // namespace jointcal
155 } // namespace lsst
156 
157 #endif // LSST_JOINTCAL_SIMPLE_PHOTOMETRY_MODEL_H
Objects used as position anchors, typically USNO stars.
Definition: RefStar.h:39
virtual double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const =0
Return the on-sky transformed flux uncertainty for measuredStar on ccdImage.
#define LOG_LOGGER
Definition: Log.h:703
void offsetFittedStar(FittedStar &fittedStar, double delta) const override
Offset the appropriate flux or magnitude (by -delta).
virtual double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const =0
Compute the residual between the model applied to a star and its associated fittedStar.
Relates transform(s) to their position in the fitting matrix and allows interaction with the transfor...
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)...
virtual std::shared_ptr< afw::image::PhotoCalib > toPhotoCalib(CcdImage const &ccdImage) const =0
Return the mapping of ccdImage represented as a PhotoCalib.
void getMappingIndices(CcdImage const &ccdImage, IndexVector &indices) const override
Get how this set of parameters (of length Npar()) map into the "grand" fit.
PhotometryMappingBase * findMapping(CcdImage const &ccdImage) const override
Return the mapping associated with this ccdImage.
void computeParameterDerivatives(MeasuredStar const &measuredStar, CcdImage const &ccdImage, Eigen::VectorXd &derivatives) const override
Compute the parametric derivatives of this model.
SimplePhotometryModel(CcdImageList const &ccdImageList, LOG_LOGGER log, double errorPedestal_=0)
double getFluxErr() const
Definition: BaseStar.h:100
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 getRefError(RefStar const &refStar) const override
Return the refStar error appropriate for this model (e.g. fluxErr or magErr).
STL class.
std::size_t getTotalParameters() const override
Return the total number of parameters in this model.
A base class for image defects.
double getRefError(RefStar const &refStar) const override
Return the refStar error appropriate for this model (e.g. fluxErr or magErr).
std::unordered_map< CcdImageKey, std::unique_ptr< PhotometryMapping > > MapType
double getMag() const
Definition: BaseStar.h:103
objects measured on actual images.
Definition: MeasuredStar.h:46
double getFlux() const
Definition: BaseStar.h:96
T size(T... args)
void offsetFittedStar(FittedStar &fittedStar, double delta) const override
Offset the appropriate flux or magnitude (by -delta).
double getMagErr() const
Definition: BaseStar.h:106
STL class.
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 offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the provided amounts (by -delta).
void dump(std::ostream &stream=std::cout) const override
Dump the contents of the transforms, for debugging.
Handler of an actual image from a single CCD.
Definition: CcdImage.h:64
STL class.
SimplePhotometryModel & operator=(SimplePhotometryModel const &)=delete
Photometric response model which has a single photometric factor per CcdImage.
The objects which have been measured several times.
Definition: FittedStar.h:61
virtual double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const =0
Return the on-sky transformed flux for measuredStar on ccdImage.
void freezeErrorTransform() override
Once this routine has been called, the error transform is not modified by offsetParams().
T reserve(T... args)