LSST Applications g063fba187b+66a50001ff,g0f08755f38+1a22dc2551,g1653933729+a8ce1bb630,g168dd56ebc+a8ce1bb630,g217e2c1bcf+12e87a5bd8,g246886dfd9+466c7b6c06,g28da252d5a+858b171e04,g2bbee38e9b+404b60ec9b,g2bc492864f+404b60ec9b,g3156d2b45e+6e55a43351,g347aa1857d+404b60ec9b,g35bb328faa+a8ce1bb630,g3a166c0a6a+404b60ec9b,g3e281a1b8c+c5dd892a6c,g414038480c+6b9177ef31,g41af890bb2+70bea58702,g599934f4f4+b8c5400ca5,g781aacb6e4+a8ce1bb630,g7af13505b9+b5b9cefdb8,g80478fca09+c2997882f3,g82479be7b0+8974e6af0f,g858d7b2824+1a22dc2551,g89c8672015+f4add4ffd5,g8f1c07a47a+de51c9b0a5,g9125e01d80+a8ce1bb630,ga5288a1d22+b66f8cf76b,gb58c049af0+d64f4d3760,gc28159a63d+404b60ec9b,gcab2d0539d+66cf1de5d4,gcf0d15dbbd+12cb7e2563,gda6a2b7d83+12cb7e2563,gdaeeff99f8+1711a396fd,ge79ae78c31+404b60ec9b,gef2f8181fd+414189b318,gf0baf85859+c1f95f4921,gf0c06eb49c+1a22dc2551,gfa517265be+1a22dc2551,gfa999e8aa5+17cd334064,v28.0.0.rc2
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
lsst::jointcal::PhotometryModel Class Referenceabstract

#include <PhotometryModel.h>

Inheritance diagram for lsst::jointcal::PhotometryModel:
lsst::jointcal::ConstrainedPhotometryModel lsst::jointcal::SimplePhotometryModel lsst::jointcal::ConstrainedFluxModel lsst::jointcal::ConstrainedMagnitudeModel lsst::jointcal::SimpleFluxModel lsst::jointcal::SimpleMagnitudeModel

Public Member Functions

 PhotometryModel (LOG_LOGGER log, double errorPedestal=0)
 
virtual Eigen::Index assignIndices (std::string const &whatToFit, Eigen::Index firstIndex)=0
 Assign indices in the full matrix to the parameters being fit in the mappings, starting at firstIndex.
 
virtual void offsetParams (Eigen::VectorXd const &delta)=0
 Offset the parameters by the provided amounts (by -delta).
 
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 void freezeErrorTransform ()=0
 Once this routine has been called, the error transform is not modified by offsetParams().
 
virtual void getMappingIndices (CcdImage const &ccdImage, IndexVector &indices) const =0
 Get how this set of parameters (of length Npar()) map into the "grand" fit.
 
virtual void computeParameterDerivatives (MeasuredStar const &measuredStar, CcdImage const &ccdImage, Eigen::VectorXd &derivatives) const =0
 Compute the parametric derivatives of this model.
 
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.
 
virtual std::size_t getTotalParameters () const =0
 Return the total number of parameters in this model.
 
virtual void print (std::ostream &out) const =0
 Print a string representation of the contents of this mapping, for debugging.
 
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.
 
virtual ~PhotometryModel ()=default
 

Protected Member Functions

virtual PhotometryMappingBasefindMapping (CcdImage const &ccdImage) const =0
 Return a pointer to the mapping associated with this ccdImage.
 

Protected Attributes

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

Friends

std::ostreamoperator<< (std::ostream &s, PhotometryModel const &model)
 

Detailed Description

Definition at line 45 of file PhotometryModel.h.

Constructor & Destructor Documentation

◆ PhotometryModel()

lsst::jointcal::PhotometryModel::PhotometryModel ( LOG_LOGGER log,
double errorPedestal = 0 )
inline
Parameters
logLogger to send messages to, to keep names consistent when logging.
errorPedestalPedestal on flux/magnitude error (percent of flux or delta magnitude).

Definition at line 51 of file PhotometryModel.h.

LOG_LOGGER _log
lsst.logging instance, to be created by a subclass so that messages have consistent name.
T move(T... args)

◆ ~PhotometryModel()

virtual lsst::jointcal::PhotometryModel::~PhotometryModel ( )
virtualdefault

Member Function Documentation

◆ assignIndices()

virtual Eigen::Index lsst::jointcal::PhotometryModel::assignIndices ( std::string const & whatToFit,
Eigen::Index firstIndex )
pure virtual

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.

Implemented in lsst::jointcal::ConstrainedPhotometryModel, and lsst::jointcal::SimplePhotometryModel.

◆ checkPositiveOnBBox()

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

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
#define LOGLS_WARN(logger, message)
Log a warn-level message using an iostream-based interface.
Definition Log.h:659
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()

virtual void lsst::jointcal::PhotometryModel::computeParameterDerivatives ( MeasuredStar const & measuredStar,
CcdImage const & ccdImage,
Eigen::VectorXd & derivatives ) const
pure virtual

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.

Implemented in lsst::jointcal::ConstrainedPhotometryModel, and lsst::jointcal::SimplePhotometryModel.

◆ computeRefResidual()

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

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 virtual

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()

virtual PhotometryMappingBase * lsst::jointcal::PhotometryModel::findMapping ( CcdImage const & ccdImage) const
protectedpure virtual

Return a pointer to the mapping associated with this ccdImage.

Implemented in lsst::jointcal::ConstrainedPhotometryModel, and lsst::jointcal::SimplePhotometryModel.

◆ freezeErrorTransform()

virtual void lsst::jointcal::PhotometryModel::freezeErrorTransform ( )
pure virtual

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.

Implemented in lsst::jointcal::ConstrainedPhotometryModel, and lsst::jointcal::SimplePhotometryModel.

◆ getErrorPedestal()

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

Definition at line 197 of file PhotometryModel.h.

197{ return errorPedestal; }

◆ getMapping()

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

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()

virtual void lsst::jointcal::PhotometryModel::getMappingIndices ( CcdImage const & ccdImage,
IndexVector & indices ) const
pure virtual

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.

Implemented in lsst::jointcal::ConstrainedPhotometryModel, and lsst::jointcal::SimplePhotometryModel.

◆ getNpar()

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

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 virtual

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()

virtual std::size_t lsst::jointcal::PhotometryModel::getTotalParameters ( ) const
pure virtual

Return the total number of parameters in this model.

Implemented in lsst::jointcal::ConstrainedPhotometryModel, and lsst::jointcal::SimplePhotometryModel.

◆ offsetFittedStar()

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

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()

virtual void lsst::jointcal::PhotometryModel::offsetParams ( Eigen::VectorXd const & delta)
pure virtual

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

Implemented in lsst::jointcal::ConstrainedPhotometryModel, and lsst::jointcal::SimplePhotometryModel.

◆ print()

virtual void lsst::jointcal::PhotometryModel::print ( std::ostream & out) const
pure virtual

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.

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

◆ toPhotoCalib()

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

◆ transform()

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

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 virtual

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
inline

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
inline

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

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.

Friends And Related Symbol Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream & s,
PhotometryModel const & model )
friend

Definition at line 192 of file PhotometryModel.h.

192 {
193 model.print(s);
194 return s;
195 }
model(data, psfmodels, sources)

Member Data Documentation

◆ _log

LOG_LOGGER lsst::jointcal::PhotometryModel::_log
protected

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

Definition at line 224 of file PhotometryModel.h.

◆ errorPedestal

double lsst::jointcal::PhotometryModel::errorPedestal
protected

Definition at line 227 of file PhotometryModel.h.


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