LSST Applications  21.0.0-172-gfb10e10a+18fedfabac,22.0.0+297cba6710,22.0.0+80564b0ff1,22.0.0+8d77f4f51a,22.0.0+a28f4c53b1,22.0.0+dcf3732eb2,22.0.1-1-g7d6de66+2a20fdde0d,22.0.1-1-g8e32f31+297cba6710,22.0.1-1-geca5380+7fa3b7d9b6,22.0.1-12-g44dc1dc+2a20fdde0d,22.0.1-15-g6a90155+515f58c32b,22.0.1-16-g9282f48+790f5f2caa,22.0.1-2-g92698f7+dcf3732eb2,22.0.1-2-ga9b0f51+7fa3b7d9b6,22.0.1-2-gd1925c9+bf4f0e694f,22.0.1-24-g1ad7a390+a9625a72a8,22.0.1-25-g5bf6245+3ad8ecd50b,22.0.1-25-gb120d7b+8b5510f75f,22.0.1-27-g97737f7+2a20fdde0d,22.0.1-32-gf62ce7b1+aa4237961e,22.0.1-4-g0b3f228+2a20fdde0d,22.0.1-4-g243d05b+871c1b8305,22.0.1-4-g3a563be+32dcf1063f,22.0.1-4-g44f2e3d+9e4ab0f4fa,22.0.1-42-gca6935d93+ba5e5ca3eb,22.0.1-5-g15c806e+85460ae5f3,22.0.1-5-g58711c4+611d128589,22.0.1-5-g75bb458+99c117b92f,22.0.1-6-g1c63a23+7fa3b7d9b6,22.0.1-6-g50866e6+84ff5a128b,22.0.1-6-g8d3140d+720564cf76,22.0.1-6-gd805d02+cc5644f571,22.0.1-8-ge5750ce+85460ae5f3,master-g6e05de7fdc+babf819c66,master-g99da0e417a+8d77f4f51a,w.2021.48
LSST Data Management Base Package
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
lsst::jointcal::SimplePhotometryModel Class Referenceabstract

Photometric response model which has a single photometric factor per CcdImage. More...

#include <SimplePhotometryModel.h>

Inheritance diagram for lsst::jointcal::SimplePhotometryModel:
lsst::jointcal::PhotometryModel lsst::jointcal::SimpleFluxModel lsst::jointcal::SimpleMagnitudeModel

Public Member Functions

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

Protected Types

typedef std::unordered_map< CcdImageKey, std::unique_ptr< PhotometryMapping > > MapType
 

Protected Member Functions

PhotometryMappingBasefindMapping (CcdImage const &ccdImage) const override
 Return the mapping associated with this ccdImage. More...
 

Protected Attributes

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

Detailed Description

Photometric response model which has a single photometric factor per CcdImage.

Definition at line 43 of file SimplePhotometryModel.h.

Member Typedef Documentation

◆ MapType

Definition at line 81 of file SimplePhotometryModel.h.

Constructor & Destructor Documentation

◆ SimplePhotometryModel() [1/3]

lsst::jointcal::SimplePhotometryModel::SimplePhotometryModel ( CcdImageList const &  ccdImageList,
LOG_LOGGER  log,
double  errorPedestal = 0 
)
inline

Definition at line 45 of file SimplePhotometryModel.h.

47  _myMap.reserve(ccdImageList.size());
48  }
PhotometryModel(LOG_LOGGER log, double errorPedestal=0)
T reserve(T... args)

◆ SimplePhotometryModel() [2/3]

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

No copy or move: there is only ever one instance of a given model.

◆ SimplePhotometryModel() [3/3]

lsst::jointcal::SimplePhotometryModel::SimplePhotometryModel ( SimplePhotometryModel &&  )
delete

◆ ~SimplePhotometryModel()

lsst::jointcal::SimplePhotometryModel::~SimplePhotometryModel ( )
default

Member Function Documentation

◆ assignIndices()

Eigen::Index lsst::jointcal::SimplePhotometryModel::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 38 of file SimplePhotometryModel.cc.

38  {
39  Eigen::Index ipar = firstIndex;
40  for (auto const &i : _myMap) {
41  auto mapping = i.second.get();
42  mapping->setIndex(ipar);
43  ipar += mapping->getNpar();
44  }
45  return ipar;
46 }

◆ 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
double x
#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
The base class for handling stars. Used by all matching routines.
Definition: BaseStar.h:51
Sources measured on images.
Definition: MeasuredStar.h:51
LOG_LOGGER _log
lsst.logging instance, to be created by a subclass so that messages have consistent name.
virtual double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const =0
Return the on-sky transformed flux for measuredStar on ccdImage.

◆ computeParameterDerivatives()

void lsst::jointcal::SimplePhotometryModel::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 76 of file SimplePhotometryModel.cc.

78  {
79  auto mapping = findMapping(ccdImage);
80  mapping->computeParameterDerivatives(measuredStar, measuredStar.getInstFlux(), derivatives);
81 }
PhotometryMappingBase * findMapping(CcdImage const &ccdImage) const override
Return 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::SimpleMagnitudeModel, lsst::jointcal::SimpleFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, and lsst::jointcal::ConstrainedFluxModel.

◆ 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::SimpleMagnitudeModel, lsst::jointcal::SimpleFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, and lsst::jointcal::ConstrainedFluxModel.

◆ findMapping()

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

Return the mapping associated with this ccdImage.

Implements lsst::jointcal::PhotometryModel.

Definition at line 91 of file SimplePhotometryModel.cc.

91  {
92  auto i = _myMap.find(ccdImage.getHashKey());
93  if (i == _myMap.end())
94  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
95  "SimplePhotometryModel cannot find CcdImage " + ccdImage.getName());
96  return i->second.get();
97 }
#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::SimplePhotometryModel::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 55 of file SimplePhotometryModel.cc.

55  {
56  for (auto &i : _myMap) {
57  i.second->freezeErrorTransform();
58  }
59 }

◆ getErrorPedestal()

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

Definition at line 196 of file PhotometryModel.h.

196 { return errorPedestal; }

◆ getMapping()

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

Get the mapping associated with ccdImage.

Definition at line 156 of file PhotometryModel.h.

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

◆ getMappingIndices()

void lsst::jointcal::SimplePhotometryModel::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 61 of file SimplePhotometryModel.cc.

62  {
63  auto mapping = findMapping(ccdImage);
64  if (indices.size() < mapping->getNpar()) indices.resize(mapping->getNpar());
65  indices[0] = mapping->getIndex();
66 }
T resize(T... args)
T size(T... args)

◆ 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 153 of file PhotometryModel.h.

153 { 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::SimpleMagnitudeModel, lsst::jointcal::SimpleFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, and lsst::jointcal::ConstrainedFluxModel.

◆ getTotalParameters()

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

Return the total number of parameters in this model.

Implements lsst::jointcal::PhotometryModel.

Definition at line 68 of file SimplePhotometryModel.cc.

68  {
69  std::size_t total = 0;
70  for (auto &i : _myMap) {
71  total += i.second->getNpar();
72  }
73  return total;
74 }

◆ 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::SimpleMagnitudeModel, lsst::jointcal::SimpleFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, and lsst::jointcal::ConstrainedFluxModel.

◆ offsetParams()

void lsst::jointcal::SimplePhotometryModel::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 48 of file SimplePhotometryModel.cc.

48  {
49  for (auto &i : _myMap) {
50  auto mapping = i.second.get();
51  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
52  }
53 }

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ print()

void lsst::jointcal::SimplePhotometryModel::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.

Reimplemented in lsst::jointcal::SimpleMagnitudeModel, and lsst::jointcal::SimpleFluxModel.

Definition at line 83 of file SimplePhotometryModel.cc.

83  {
84  for (auto &i : _myMap) {
85  out << i.first << ": ";
86  i.second->print(out);
87  out << std::endl;
88  }
89 }
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::SimpleMagnitudeModel, lsst::jointcal::SimpleFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, and lsst::jointcal::ConstrainedFluxModel.

◆ 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::SimpleMagnitudeModel, lsst::jointcal::SimpleFluxModel, lsst::jointcal::ConstrainedMagnitudeModel, and lsst::jointcal::ConstrainedFluxModel.

◆ 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 199 of file PhotometryModel.h.

199  {
200  if (errorPedestal == 0) {
201  return measuredStar.getInstFluxErr();
202  } else {
203  return std::hypot(measuredStar.getInstFluxErr(), measuredStar.getInstFlux() * errorPedestal);
204  }
205  }
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 208 of file PhotometryModel.h.

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

◆ 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

◆ _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 223 of file PhotometryModel.h.

◆ _myMap

MapType lsst::jointcal::SimplePhotometryModel::_myMap
protected

Definition at line 82 of file SimplePhotometryModel.h.

◆ errorPedestal

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

Definition at line 226 of file PhotometryModel.h.


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