LSST Applications g02d81e74bb+86cf3d8bc9,g180d380827+7a4e862ed4,g2079a07aa2+86d27d4dc4,g2305ad1205+e1ca1c66fa,g29320951ab+012e1474a1,g295015adf3+341ea1ce94,g2bbee38e9b+0e5473021a,g337abbeb29+0e5473021a,g33d1c0ed96+0e5473021a,g3a166c0a6a+0e5473021a,g3ddfee87b4+c429d67c83,g48712c4677+f88676dd22,g487adcacf7+27e1e21933,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+b41db86c35,g5a732f18d5+53520f316c,g64a986408d+86cf3d8bc9,g858d7b2824+86cf3d8bc9,g8a8a8dda67+585e252eca,g99cad8db69+84912a7fdc,g9ddcbc5298+9a081db1e4,ga1e77700b3+15fc3df1f7,ga8c6da7877+a2b54eae19,gb0e22166c9+60f28cb32d,gba4ed39666+c2a2e4ac27,gbb8dafda3b+6681f309db,gc120e1dc64+f0fcc2f6d8,gc28159a63d+0e5473021a,gcf0d15dbbd+c429d67c83,gdaeeff99f8+f9a426f77a,ge6526c86ff+0433e6603d,ge79ae78c31+0e5473021a,gee10cc3b42+585e252eca,gff1a9f87cc+86cf3d8bc9,w.2024.17
LSST Data Management Base Package
Loading...
Searching...
No Matches
Classes | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
lsst::jointcal::ConstrainedPhotometryModel Class Referenceabstract

Photometry model with constraints, \(M(x,y) = M_CCD(x,y)*M_visit(u,v)\). More...

#include <ConstrainedPhotometryModel.h>

Inheritance diagram for lsst::jointcal::ConstrainedPhotometryModel:
lsst::jointcal::PhotometryModel lsst::jointcal::ConstrainedFluxModel lsst::jointcal::ConstrainedMagnitudeModel

Classes

struct  PrepPhotoCalib
 To hold the return of prepPhotoCalib. More...
 

Public Member Functions

 ConstrainedPhotometryModel (CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox, LOG_LOGGER log, int visitOrder=7, double errorPedestal=0)
 Construct a constrained photometry model.
 
 ConstrainedPhotometryModel (ConstrainedPhotometryModel const &)=delete
 No copy or move: there is only ever one instance of a given model (i.e. per ccd+visit)
 
 ConstrainedPhotometryModel (ConstrainedPhotometryModel &&)=delete
 
ConstrainedPhotometryModeloperator= (ConstrainedPhotometryModel const &)=delete
 
ConstrainedPhotometryModeloperator= (ConstrainedPhotometryModel &&)=delete
 
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 freezeErrorTransform () override
 Once this routine has been called, the error transform is not modified by offsetParams().
 
void getMappingIndices (CcdImage const &ccdImage, IndexVector &indices) const override
 Get how this set of parameters (of length Npar()) map into the "grand" fit.
 
std::size_t getTotalParameters () const override
 Return the total number of parameters in this model.
 
void computeParameterDerivatives (MeasuredStar const &measuredStar, CcdImage const &ccdImage, Eigen::VectorXd &derivatives) const override
 Compute the parametric derivatives of this model.
 
void print (std::ostream &out) const override
 Print a string representation of the contents of this mapping, for debugging.
 
 ~ConstrainedPhotometryModel ()=default
 
virtual void offsetFittedStar (FittedStar &fittedStar, double delta) const =0
 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.
 
virtual double transform (CcdImage const &ccdImage, MeasuredStar const &measuredStar) const =0
 Return the on-sky transformed flux for measuredStar on ccdImage.
 
virtual double transformError (CcdImage const &ccdImage, MeasuredStar const &measuredStar) const =0
 Return the on-sky transformed flux uncertainty for measuredStar on ccdImage.
 
virtual double getRefError (RefStar const &refStar) const =0
 Return the refStar error appropriate for this model (e.g. fluxErr or magErr).
 
virtual double computeRefResidual (FittedStar const &fittedStar, RefStar const &refStar) const =0
 Return the fittedStar - refStar residual appropriate for this model (e.g. flux - flux or mag - mag).
 
virtual std::shared_ptr< afw::image::PhotoCalibtoPhotoCalib (CcdImage const &ccdImage) const =0
 Return the mapping of ccdImage represented as a PhotoCalib.
 
std::size_t getNpar (CcdImage const &ccdImage) const
 Return the number of parameters in the mapping of CcdImage.
 
PhotometryMappingBase const & getMapping (CcdImage const &ccdImage) const
 Get the mapping associated with ccdImage.
 
bool validate (CcdImageList const &ccdImageList, int ndof) const
 Return true if this is a "reasonable" model.
 
bool checkPositiveOnBBox (CcdImage const &ccdImage) const
 Check that the model is positive on the ccdImage bbox.
 
double getErrorPedestal () const
 
double tweakFluxError (jointcal::MeasuredStar const &measuredStar) const
 Add a fraction of the instrumental flux to the instrumental flux error, in quadrature.
 
double tweakMagnitudeError (jointcal::MeasuredStar const &measuredStar) const
 Add a small magnitude offset to the "instrumental magnitude" error, in quadrature.
 

Protected Types

using MapType = std::unordered_map<CcdImageKey, std::unique_ptr<ChipVisitPhotometryMapping>>
 
using VisitMapType = std::map<VisitIdType, std::shared_ptr<PhotometryMapping>>
 
using ChipMapType = std::map<CcdIdType, std::shared_ptr<PhotometryMapping>>
 

Protected Member Functions

PhotometryMappingBasefindMapping (CcdImage const &ccdImage) const override
 Return a pointer to the mapping associated with this ccdImage.
 
template<class ChipTransform , class VisitTransform , class ChipVisitMapping >
void initialize (CcdImageList const &ccdImageList, geom::Box2D const &focalPlaneBBox, int visitOrder)
 Initialize the chip, visit, and chipVisit mappings by creating appropriate transforms and mappings.
 
virtual double initialChipCalibration (std::shared_ptr< afw::image::PhotoCalib const > photoCalib)=0
 Return the initial calibration to use from this photoCalib.
 
PrepPhotoCalib prepPhotoCalib (CcdImage const &ccdImage) const
 Helper for preparing toPhotoCalib()
 

Protected Attributes

MapType _chipVisitMap
 
VisitMapType _visitMap
 
ChipMapType _chipMap
 
LOG_LOGGER _log
 lsst.logging instance, to be created by a subclass so that messages have consistent name.
 
double errorPedestal
 

Detailed Description

Photometry model with constraints, \(M(x,y) = M_CCD(x,y)*M_visit(u,v)\).

This model consists of the following components:

Because this model's parameters are degenerate under multiplication by a constant, \(M=(a*M_CCD)*(1/a*M_visit)\), we hold one CCD's zero point fixed to remove that degeneracy.

Definition at line 49 of file ConstrainedPhotometryModel.h.

Member Typedef Documentation

◆ ChipMapType

Definition at line 111 of file ConstrainedPhotometryModel.h.

◆ MapType

Definition at line 104 of file ConstrainedPhotometryModel.h.

◆ VisitMapType

Definition at line 108 of file ConstrainedPhotometryModel.h.

Constructor & Destructor Documentation

◆ ConstrainedPhotometryModel() [1/3]

lsst::jointcal::ConstrainedPhotometryModel::ConstrainedPhotometryModel ( CcdImageList const & ccdImageList,
geom::Box2D const & focalPlaneBBox,
LOG_LOGGER log,
int visitOrder = 7,
double errorPedestal = 0 )
inlineexplicit

Construct a constrained photometry model.

Parameters
ccdImageListThe list of CCDImages to construct the model for.
focalPlaneBBoxThe bounding box of the camera's focal plane, defining the domain of the visit polynomial.
[in]logAn lsst::log::Log instance to log messages to.
[in]visitOrderThe order of the visit polynomial.
[in]errorPedestalA pedestal in flux or magnitude to apply to all MeasuredStar flux errors.

Definition at line 61 of file ConstrainedPhotometryModel.h.

63 : PhotometryModel(std::move(log), errorPedestal), _fittingChips(false), _fittingVisits(false) {
64 _chipVisitMap.reserve(ccdImageList.size());
65 }
PhotometryModel(LOG_LOGGER log, double errorPedestal=0)
T move(T... args)
T reserve(T... args)

◆ ConstrainedPhotometryModel() [2/3]

lsst::jointcal::ConstrainedPhotometryModel::ConstrainedPhotometryModel ( ConstrainedPhotometryModel const & )
delete

No copy or move: there is only ever one instance of a given model (i.e. per ccd+visit)

◆ ConstrainedPhotometryModel() [3/3]

lsst::jointcal::ConstrainedPhotometryModel::ConstrainedPhotometryModel ( ConstrainedPhotometryModel && )
delete

◆ ~ConstrainedPhotometryModel()

lsst::jointcal::ConstrainedPhotometryModel::~ConstrainedPhotometryModel ( )
default

Member Function Documentation

◆ assignIndices()

Eigen::Index lsst::jointcal::ConstrainedPhotometryModel::assignIndices ( std::string const & whatToFit,
Eigen::Index firstIndex )
overridevirtual

Assign indices in the full matrix to the parameters being fit in the mappings, starting at firstIndex.

Parameters
[in]whatToFitString containing parameters to fit.
[in]firstIndexIndex to start assigning at.
Returns
The highest assigned index.

Implements lsst::jointcal::PhotometryModel.

Definition at line 45 of file ConstrainedPhotometryModel.cc.

46 {
47 Eigen::Index index = firstIndex;
48 if (whatToFit.find("Model") == std::string::npos) {
49 LOGLS_WARN(_log, "assignIndices was called and Model is *not* in whatToFit");
50 return index;
51 }
52
53 // If we got here, "Model" is definitely in whatToFit.
54 _fittingChips = (whatToFit.find("ModelChip") != std::string::npos);
55 _fittingVisits = (whatToFit.find("ModelVisit") != std::string::npos);
56 // If nothing more than "Model" is specified, it means fit everything.
57 if ((!_fittingChips) && (!_fittingVisits)) {
58 _fittingChips = _fittingVisits = true;
59 }
60
61 if (_fittingChips) {
62 for (auto &idMapping : _chipMap) {
63 auto mapping = idMapping.second.get();
64 // Don't assign indices for fixed parameters.
65 if (mapping->isFixed()) continue;
66 mapping->setIndex(index);
67 index += mapping->getNpar();
68 }
69 }
70 if (_fittingVisits) {
71 for (auto &idMapping : _visitMap) {
72 auto mapping = idMapping.second.get();
73 mapping->setIndex(index);
74 index += mapping->getNpar();
75 }
76 }
77 for (auto &idMapping : _chipVisitMap) {
78 idMapping.second->setWhatToFit(_fittingChips, _fittingVisits);
79 }
80 return index;
81}
#define LOGLS_WARN(logger, message)
Log a warn-level message using an iostream-based interface.
Definition Log.h:659
LOG_LOGGER _log
lsst.logging instance, to be created by a subclass so that messages have consistent name.

◆ checkPositiveOnBBox()

bool lsst::jointcal::PhotometryModel::checkPositiveOnBBox ( CcdImage const & ccdImage) const
inherited

Check that the model is positive on the ccdImage bbox.

Parameters
ccdImageThe ccdImage to test.
Returns
True if the image is positive on a sampling of points of the ccdImage bbox.

Definition at line 49 of file PhotometryModel.cc.

49 {
50 bool check = true;
51 auto bbox = ccdImage.getImageFrame();
52 for (auto const &x : {bbox.xMin, bbox.getCenter().x, bbox.xMax})
53 for (auto const &y : {bbox.yMin, bbox.getCenter().y, bbox.yMax}) {
54 // flux, fluxErr, instFlux, instFluxErr all 1
56 star.setInstFluxAndErr(1, 1);
57 double result = transform(ccdImage, star);
58 // Don't short circuit so that we log every place the model is negative.
59 if (result < 0) {
60 LOGLS_WARN(_log, "CcdImage " << ccdImage.getName() << " is negative at (" << x << "," << y
61 << "): " << result);
62 check = false;
63 }
64 }
65 return check;
66}
py::object result
Definition _schema.cc:429
AmpInfoBoxKey bbox
Definition Amplifier.cc:117
int y
Definition SpanSet.cc:48
table::Key< int > transform
The base class for handling stars. Used by all matching routines.
Definition BaseStar.h:51
Sources measured on images.

◆ computeParameterDerivatives()

void lsst::jointcal::ConstrainedPhotometryModel::computeParameterDerivatives ( MeasuredStar const & measuredStar,
CcdImage const & ccdImage,
Eigen::VectorXd & derivatives ) const
overridevirtual

Compute the parametric derivatives of this model.

Parameters
[in]measuredStarThe measured star with the position and flux to compute at.
[in]ccdImageThe ccdImage containing the measured star, to find the correct mapping.
[out]derivativesThe computed derivatives. Must be pre-allocated to the correct size.

Implements lsst::jointcal::PhotometryModel.

Definition at line 125 of file ConstrainedPhotometryModel.cc.

127 {
128 auto mapping = findMapping(ccdImage);
129 mapping->computeParameterDerivatives(measuredStar, measuredStar.getInstFlux(), derivatives);
130}
PhotometryMappingBase * findMapping(CcdImage const &ccdImage) const override
Return a pointer to the mapping associated with this ccdImage.

◆ computeRefResidual()

virtual double lsst::jointcal::PhotometryModel::computeRefResidual ( FittedStar const & fittedStar,
RefStar const & refStar ) const
pure virtualinherited

Return the fittedStar - refStar residual appropriate for this model (e.g. flux - flux or mag - mag).

Implemented in lsst::jointcal::ConstrainedFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, lsst::jointcal::SimpleFluxModel, and lsst::jointcal::SimpleMagnitudeModel.

◆ computeResidual()

virtual double lsst::jointcal::PhotometryModel::computeResidual ( CcdImage const & ccdImage,
MeasuredStar const & measuredStar ) const
pure virtualinherited

Compute the residual between the model applied to a star and its associated fittedStar.

\[ residual = Model(measuredStar) - fittedStar \]

Parameters
ccdImageThe ccdImage where measuredStar resides.
measuredStarThe measured star position to compute the residual of.
Returns
The residual.

Implemented in lsst::jointcal::ConstrainedFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, lsst::jointcal::SimpleFluxModel, and lsst::jointcal::SimpleMagnitudeModel.

◆ findMapping()

PhotometryMappingBase * lsst::jointcal::ConstrainedPhotometryModel::findMapping ( CcdImage const & ccdImage) const
overrideprotectedvirtual

Return a pointer to the mapping associated with this ccdImage.

Implements lsst::jointcal::PhotometryModel.

Definition at line 168 of file ConstrainedPhotometryModel.cc.

168 {
169 auto idMapping = _chipVisitMap.find(ccdImage.getHashKey());
170 if (idMapping == _chipVisitMap.end())
171 throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
172 "ConstrainedPhotometryModel cannot find CcdImage " + ccdImage.getName());
173 return idMapping->second.get();
174}
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition Exception.h:48
T end(T... args)
T find(T... args)

◆ freezeErrorTransform()

void lsst::jointcal::ConstrainedPhotometryModel::freezeErrorTransform ( )
overridevirtual

Once this routine has been called, the error transform is not modified by offsetParams().

The routine can be called when the mappings are roughly in place. After the call, the transformations used to propagate errors are no longer affected when updating the mappings. This allows an exactly linear fit, which can be necessary for some model+data combinations.

Implements lsst::jointcal::PhotometryModel.

Definition at line 100 of file ConstrainedPhotometryModel.cc.

100 {
101 for (auto &idMapping : _chipMap) {
102 idMapping.second->freezeErrorTransform();
103 }
104 for (auto &idMapping : _visitMap) {
105 idMapping.second->freezeErrorTransform();
106 }
107}

◆ getErrorPedestal()

double lsst::jointcal::PhotometryModel::getErrorPedestal ( ) const
inlineinherited

Definition at line 197 of file PhotometryModel.h.

197{ return errorPedestal; }

◆ getMapping()

PhotometryMappingBase const & lsst::jointcal::PhotometryModel::getMapping ( CcdImage const & ccdImage) const
inlineinherited

Get the mapping associated with ccdImage.

Definition at line 157 of file PhotometryModel.h.

157 {
158 return *(findMapping(ccdImage));
159 }
virtual PhotometryMappingBase * findMapping(CcdImage const &ccdImage) const =0
Return a pointer to the mapping associated with this ccdImage.

◆ getMappingIndices()

void lsst::jointcal::ConstrainedPhotometryModel::getMappingIndices ( CcdImage const & ccdImage,
IndexVector & indices ) const
overridevirtual

Get how this set of parameters (of length Npar()) map into the "grand" fit.

Parameters
[in]ccdImageThe ccdImage to look up.
[out]indicesThe indices of the mapping associated with ccdImage.

Implements lsst::jointcal::PhotometryModel.

Definition at line 109 of file ConstrainedPhotometryModel.cc.

109 {
110 auto mapping = findMapping(ccdImage);
111 mapping->getMappingIndices(indices);
112}

◆ getNpar()

std::size_t lsst::jointcal::PhotometryModel::getNpar ( CcdImage const & ccdImage) const
inlineinherited

Return the number of parameters in the mapping of CcdImage.

Definition at line 154 of file PhotometryModel.h.

154{ return findMapping(ccdImage)->getNpar(); }
virtual std::size_t getNpar() const =0
Number of total parameters in this mapping.

◆ getRefError()

virtual double lsst::jointcal::PhotometryModel::getRefError ( RefStar const & refStar) const
pure virtualinherited

Return the refStar error appropriate for this model (e.g. fluxErr or magErr).

Implemented in lsst::jointcal::ConstrainedFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, lsst::jointcal::SimpleFluxModel, and lsst::jointcal::SimpleMagnitudeModel.

◆ getTotalParameters()

std::size_t lsst::jointcal::ConstrainedPhotometryModel::getTotalParameters ( ) const
overridevirtual

Return the total number of parameters in this model.

Implements lsst::jointcal::PhotometryModel.

Definition at line 114 of file ConstrainedPhotometryModel.cc.

114 {
115 std::size_t total = 0;
116 for (auto &idMapping : _chipMap) {
117 total += idMapping.second->getNpar();
118 }
119 for (auto &idMapping : _visitMap) {
120 total += idMapping.second->getNpar();
121 }
122 return total;
123}

◆ initialChipCalibration()

virtual double lsst::jointcal::ConstrainedPhotometryModel::initialChipCalibration ( std::shared_ptr< afw::image::PhotoCalib const > photoCalib)
protectedpure virtual

Return the initial calibration to use from this photoCalib.

Implemented in lsst::jointcal::ConstrainedFluxModel, and lsst::jointcal::ConstrainedMagnitudeModel.

◆ initialize()

template<class ChipTransform , class VisitTransform , class ChipVisitMapping >
template void lsst::jointcal::ConstrainedPhotometryModel::initialize< MagnitudeTransformSpatiallyInvariant, MagnitudeTransformChebyshev, ChipVisitMagnitudeMapping > ( CcdImageList const & ccdImageList,
geom::Box2D const & focalPlaneBBox,
int visitOrder )
protected

Initialize the chip, visit, and chipVisit mappings by creating appropriate transforms and mappings.

Definition at line 177 of file ConstrainedPhotometryModel.cc.

178 {
179 // keep track of which chip we want to constrain (the one closest to the middle of the focal plane)
180 double minRadius2 = std::numeric_limits<double>::infinity();
181 CcdIdType constrainedChip = -1;
182
183 // First initialize all visit and ccd transforms, before we make the ccdImage mappings.
184 for (auto const &ccdImage : ccdImageList) {
185 auto visit = ccdImage->getVisit();
186 auto chip = ccdImage->getCcdId();
187 auto visitPair = _visitMap.find(visit);
188 auto chipPair = _chipMap.find(chip);
189
190 // If the chip is not in the map, add it, otherwise continue.
191 if (chipPair == _chipMap.end()) {
192 auto center = ccdImage->getDetector()->getCenter(afw::cameraGeom::FOCAL_PLANE);
193 double radius2 = std::pow(center.getX(), 2) + std::pow(center.getY(), 2);
194 if (radius2 < minRadius2) {
195 minRadius2 = radius2;
196 constrainedChip = chip;
197 }
198 auto photoCalib = ccdImage->getPhotoCalib();
199 // Use the single-frame processing calibration from the PhotoCalib as the default.
200 auto chipTransform = std::make_unique<ChipTransform>(initialChipCalibration(photoCalib));
201 _chipMap[chip] = std::make_shared<PhotometryMapping>(std::move(chipTransform));
202 }
203 // If the visit is not in the map, add it, otherwise continue.
204 if (visitPair == _visitMap.end()) {
205 auto visitTransform = std::make_unique<VisitTransform>(visitOrder, focalPlaneBBox);
206 _visitMap[visit] = std::make_shared<PhotometryMapping>(std::move(visitTransform));
207 }
208 }
209
210 // Fix one chip mapping, to remove the degeneracy from the system.
211 _chipMap.at(constrainedChip)->setFixed(true);
212
213 // Now create the ccdImage mappings, which are combinations of the chip/visit mappings above.
214 for (auto const &ccdImage : ccdImageList) {
215 auto visit = ccdImage->getVisit();
216 auto chip = ccdImage->getCcdId();
217 _chipVisitMap.emplace(ccdImage->getHashKey(),
218 std::make_unique<ChipVisitMapping>(_chipMap[chip], _visitMap[visit]));
219 }
220 LOGLS_INFO(_log, "Got " << _chipMap.size() << " chip mappings and " << _visitMap.size()
221 << " visit mappings; holding chip " << constrainedChip << " fixed ("
222 << getTotalParameters() << " total parameters).");
223 LOGLS_DEBUG(_log, "CcdImage map has " << _chipVisitMap.size() << " mappings, with "
224 << _chipVisitMap.bucket_count() << " buckets and a load factor of "
226}
#define LOGLS_INFO(logger, message)
Log a info-level message using an iostream-based interface.
Definition Log.h:639
#define LOGLS_DEBUG(logger, message)
Log a debug-level message using an iostream-based interface.
Definition Log.h:619
T at(T... args)
T bucket_count(T... args)
std::size_t getTotalParameters() const override
Return the total number of parameters in this model.
virtual double initialChipCalibration(std::shared_ptr< afw::image::PhotoCalib const > photoCalib)=0
Return the initial calibration to use from this photoCalib.
T emplace(T... args)
T infinity(T... args)
T load_factor(T... args)
CameraSys const FOCAL_PLANE
Focal plane coordinates: Position on a 2-d planar approximation to the focal plane (x,...
Definition CameraSys.cc:30
T pow(T... args)
T size(T... args)

◆ offsetFittedStar()

virtual void lsst::jointcal::PhotometryModel::offsetFittedStar ( FittedStar & fittedStar,
double delta ) const
pure virtualinherited

Offset the appropriate flux or magnitude (by -delta).

Parameters
fittedStarThe star to update.
deltaThe amount to update by.

Implemented in lsst::jointcal::ConstrainedFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, lsst::jointcal::SimpleFluxModel, and lsst::jointcal::SimpleMagnitudeModel.

◆ offsetParams()

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

Offset the parameters by the provided amounts (by -delta).

The shifts are applied according to the indices given in assignIndices.

Parameters
[in]deltavector of offsets to apply

Implements lsst::jointcal::PhotometryModel.

Definition at line 83 of file ConstrainedPhotometryModel.cc.

83 {
84 if (_fittingChips) {
85 for (auto &idMapping : _chipMap) {
86 auto mapping = idMapping.second.get();
87 // Don't offset indices for fixed parameters.
88 if (mapping->isFixed()) continue;
89 mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
90 }
91 }
92 if (_fittingVisits) {
93 for (auto &idMapping : _visitMap) {
94 auto mapping = idMapping.second.get();
95 mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
96 }
97 }
98}

◆ operator=() [1/2]

ConstrainedPhotometryModel & lsst::jointcal::ConstrainedPhotometryModel::operator= ( ConstrainedPhotometryModel && )
delete

◆ operator=() [2/2]

ConstrainedPhotometryModel & lsst::jointcal::ConstrainedPhotometryModel::operator= ( ConstrainedPhotometryModel const & )
delete

◆ prepPhotoCalib()

ConstrainedPhotometryModel::PrepPhotoCalib lsst::jointcal::ConstrainedPhotometryModel::prepPhotoCalib ( CcdImage const & ccdImage) const
protected

Helper for preparing toPhotoCalib()

Definition at line 228 of file ConstrainedPhotometryModel.cc.

229 {
230 auto detector = ccdImage.getDetector();
231 auto ccdBBox = detector->getBBox();
232 auto *mapping = dynamic_cast<ChipVisitPhotometryMapping *>(findMapping(ccdImage));
233
234 // There should be no way in which we can get to this point and not have a ChipVisitMapping,
235 // so blow up if we don't.
236 assert(mapping != nullptr);
237 // We know it's a Chebyshev transform because we created it as such, so blow up if it's not.
238 auto visitPhotometryTransform = std::dynamic_pointer_cast<PhotometryTransformChebyshev>(
239 mapping->getVisitMapping()->getTransform());
240 assert(visitPhotometryTransform != nullptr);
241 auto focalBBox = visitPhotometryTransform->getBBox();
242
243 // Unravel our chebyshev coefficients to build an astshim::ChebyMap.
244 auto coeff_f = toChebyMapCoeffs(std::dynamic_pointer_cast<PhotometryTransformChebyshev>(
245 mapping->getVisitMapping()->getTransform()));
246 // Bounds are the bbox
247 std::vector<double> lowerBound = {focalBBox.getMinX(), focalBBox.getMinY()};
248 std::vector<double> upperBound = {focalBBox.getMaxX(), focalBBox.getMaxY()};
249 afw::geom::TransformPoint2ToGeneric visitTransform(ast::ChebyMap(coeff_f, 1, lowerBound, upperBound));
250
251 double chipConstant = mapping->getChipMapping()->getParameters()[0];
252
253 // Compute a box that covers the area of the ccd in focal plane coordinates.
254 // This is the box over which we want to compute the mean of the visit transform.
255 auto pixToFocal = detector->getTransform(afw::cameraGeom::PIXELS, afw::cameraGeom::FOCAL_PLANE);
256 geom::Box2D ccdBBoxInFocal;
257 for (auto const &point : pixToFocal->applyForward(geom::Box2D(ccdBBox).getCorners())) {
258 ccdBBoxInFocal.include(point);
259 }
260 double visitMean = visitPhotometryTransform->mean(ccdBBoxInFocal);
261
262 return {chipConstant, visitTransform, pixToFocal, visitMean};
263}
A ChebyMap is a form of Mapping which performs a Chebyshev polynomial transformation.
Definition ChebyMap.h:97
A floating-point coordinate rectangle geometry.
Definition Box.h:413
void include(Point2D const &point) noexcept
Expand this to ensure that this->contains(point).
Definition Box.cc:380
CameraSysPrefix const PIXELS
Pixel coordinates: Nominal position on the entry surface of a given detector (x, y unbinned pixels).
Definition CameraSys.cc:34
Transform< Point2Endpoint, GenericEndpoint > TransformPoint2ToGeneric
Definition Transform.h:301

◆ print()

void lsst::jointcal::ConstrainedPhotometryModel::print ( std::ostream & out) const
overridevirtual

Print a string representation of the contents of this mapping, for debugging.

This string representation can be very verbose, as it contains all of the parameters of all of the transforms in this model.

Implements lsst::jointcal::PhotometryModel.

Definition at line 154 of file ConstrainedPhotometryModel.cc.

154 {
155 for (auto &idMapping : _chipMap) {
156 out << "Sensor: " << idMapping.first << std::endl;
157 idMapping.second->print(out);
158 out << std::endl;
159 }
160 out << std::endl;
161 for (auto &idMapping : _visitMap) {
162 out << "Visit: " << idMapping.first << std::endl;
163 idMapping.second->print(out);
164 out << std::endl;
165 }
166}
T endl(T... args)

◆ toPhotoCalib()

virtual std::shared_ptr< afw::image::PhotoCalib > lsst::jointcal::PhotometryModel::toPhotoCalib ( CcdImage const & ccdImage) const
pure virtualinherited

◆ transform()

virtual double lsst::jointcal::PhotometryModel::transform ( CcdImage const & ccdImage,
MeasuredStar const & measuredStar ) const
pure virtualinherited

Return the on-sky transformed flux for measuredStar on ccdImage.

Parameters
[in]ccdImageThe ccdImage where measuredStar resides.
measuredStarThe measured star position to transform.
Returns
The on-sky flux transformed from instFlux at measuredStar's position.

Implemented in lsst::jointcal::ConstrainedFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, lsst::jointcal::SimpleFluxModel, and lsst::jointcal::SimpleMagnitudeModel.

◆ transformError()

virtual double lsst::jointcal::PhotometryModel::transformError ( CcdImage const & ccdImage,
MeasuredStar const & measuredStar ) const
pure virtualinherited

Return the on-sky transformed flux uncertainty for measuredStar on ccdImage.

Identical to transform() until freezeErrorTransform() is called.

Parameters
[in]ccdImageThe ccdImage where measuredStar resides.
measuredStarThe measured star position to transform.
Returns
The on-sky flux transformed from instFlux at measuredStar's position.

Implemented in lsst::jointcal::ConstrainedFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, lsst::jointcal::SimpleFluxModel, and lsst::jointcal::SimpleMagnitudeModel.

◆ tweakFluxError()

double lsst::jointcal::PhotometryModel::tweakFluxError ( jointcal::MeasuredStar const & measuredStar) const
inlineinherited

Add a fraction of the instrumental flux to the instrumental flux error, in quadrature.

Definition at line 200 of file PhotometryModel.h.

200 {
201 if (errorPedestal == 0) {
202 return measuredStar.getInstFluxErr();
203 } else {
204 return std::hypot(measuredStar.getInstFluxErr(), measuredStar.getInstFlux() * errorPedestal);
205 }
206 }
T hypot(T... args)

◆ tweakMagnitudeError()

double lsst::jointcal::PhotometryModel::tweakMagnitudeError ( jointcal::MeasuredStar const & measuredStar) const
inlineinherited

Add a small magnitude offset to the "instrumental magnitude" error, in quadrature.

Definition at line 209 of file PhotometryModel.h.

209 {
210 if (errorPedestal == 0) {
211 return measuredStar.getInstMagErr();
212 } else {
213 return std::hypot(measuredStar.getInstMagErr(), errorPedestal);
214 }
215 }

◆ validate()

bool lsst::jointcal::PhotometryModel::validate ( CcdImageList const & ccdImageList,
int ndof ) const
inherited

Return true if this is a "reasonable" model.

A valid photometry model is positive within each sensor's bounding box.

Parameters
ccdImageListThe ccdImages to test the model validity on.
ndofThe number of degrees of freedom in the fit, e.g. from Fitterbase.computeChi2().
Returns
True if the model is valid on all ccdImages.

Definition at line 33 of file PhotometryModel.cc.

33 {
34 bool check = true;
35 for (auto const &ccdImage : ccdImageList) {
36 // Don't short circuit so that we log every place the model is negative.
37 check &= checkPositiveOnBBox(*ccdImage);
38 }
39 if (ndof < 0) {
40 check &= false;
41 LOGLS_ERROR(_log, "This model only has "
42 << ndof << " degrees of freedom, with " << getTotalParameters()
43 << " total parameters. Reduce the model complexity (e.g. polynomial order)"
44 " to better match the number of measured sources.");
45 }
46 return check;
47}
#define LOGLS_ERROR(logger, message)
Log a error-level message using an iostream-based interface.
Definition Log.h:679
virtual std::size_t getTotalParameters() const =0
Return the total number of parameters in this model.
bool checkPositiveOnBBox(CcdImage const &ccdImage) const
Check that the model is positive on the ccdImage bbox.

Member Data Documentation

◆ _chipMap

ChipMapType lsst::jointcal::ConstrainedPhotometryModel::_chipMap
protected

Definition at line 112 of file ConstrainedPhotometryModel.h.

◆ _chipVisitMap

MapType lsst::jointcal::ConstrainedPhotometryModel::_chipVisitMap
protected

Definition at line 105 of file ConstrainedPhotometryModel.h.

◆ _log

LOG_LOGGER lsst::jointcal::PhotometryModel::_log
protectedinherited

lsst.logging instance, to be created by a subclass so that messages have consistent name.

Definition at line 224 of file PhotometryModel.h.

◆ _visitMap

VisitMapType lsst::jointcal::ConstrainedPhotometryModel::_visitMap
protected

Definition at line 109 of file ConstrainedPhotometryModel.h.

◆ errorPedestal

double lsst::jointcal::PhotometryModel::errorPedestal
protectedinherited

Definition at line 227 of file PhotometryModel.h.


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