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 | Static Public Member Functions | Protected Types | Protected Member Functions | Static Protected Member Functions | Related Functions | List of all members
lsst::afw::image::PhotoCalib Class Referencefinal

The photometric calibration of an exposure. More...

#include <PhotoCalib.h>

Inheritance diagram for lsst::afw::image::PhotoCalib:
lsst::afw::table::io::PersistableFacade< PhotoCalib > lsst::afw::typehandling::Storable lsst::afw::table::io::Persistable

Public Member Functions

 PhotoCalib (PhotoCalib const &)=default
 
 PhotoCalib (PhotoCalib &&)=default
 
PhotoCaliboperator= (PhotoCalib const &)=delete
 
PhotoCaliboperator= (PhotoCalib &&)=delete
 
 ~PhotoCalib () override=default
 
 PhotoCalib ()
 Create a empty, zeroed calibration. More...
 
 PhotoCalib (double calibrationMean, double calibrationErr=0, lsst::geom::Box2I const &bbox=lsst::geom::Box2I())
 Create a non-spatially-varying calibration. More...
 
 PhotoCalib (std::shared_ptr< afw::math::BoundedField > calibration, double calibrationErr=0)
 Create a spatially-varying calibration. More...
 
 PhotoCalib (double calibrationMean, double calibrationErr, std::shared_ptr< afw::math::BoundedField > calibration, bool isConstant)
 Create a calibration with a pre-computed mean. More...
 
double instFluxToNanojansky (double instFlux, lsst::geom::Point< double, 2 > const &point) const
 Convert instFlux in ADU to nJy at a point in the BoundedField. More...
 
double instFluxToNanojansky (double instFlux) const
 
Measurement instFluxToNanojansky (double instFlux, double instFluxErr, lsst::geom::Point< double, 2 > const &point) const
 Convert instFlux and error in instFlux (ADU) to nJy and nJy error. More...
 
Measurement instFluxToNanojansky (double instFlux, double instFluxErr) const
 
Measurement instFluxToNanojansky (const afw::table::SourceRecord &sourceRecord, std::string const &instFluxField) const
 Convert sourceRecord[instFluxField_instFlux] (ADU) at location (sourceRecord.get("x"), sourceRecord.get("y")) (pixels) to flux and flux error (in nJy). More...
 
ndarray::Array< double, 2, 2 > instFluxToNanojansky (afw::table::SourceCatalog const &sourceCatalog, std::string const &instFluxField) const
 Convert sourceCatalog[instFluxField_instFlux] (ADU) at locations (sourceCatalog.get("x"), sourceCatalog.get("y")) (pixels) to nJy. More...
 
void instFluxToNanojansky (afw::table::SourceCatalog &sourceCatalog, std::string const &instFluxField, std::string const &outField) const
 Convert sourceCatalog[instFluxField_instFlux] (ADU) at locations (sourceCatalog.get("x"), sourceCatalog.get("y")) (pixels) to nJy and write the results back to sourceCatalog[outField_mag]. More...
 
double instFluxToMagnitude (double instFlux, lsst::geom::Point< double, 2 > const &point) const
 Convert instFlux in ADU to AB magnitude. More...
 
double instFluxToMagnitude (double instFlux) const
 
Measurement instFluxToMagnitude (double instFlux, double instFluxErr, lsst::geom::Point< double, 2 > const &point) const
 Convert instFlux and error in instFlux (ADU) to AB magnitude and magnitude error. More...
 
Measurement instFluxToMagnitude (double instFlux, double instFluxErr) const
 
Measurement instFluxToMagnitude (afw::table::SourceRecord const &sourceRecord, std::string const &instFluxField) const
 Convert sourceRecord[instFluxField_instFlux] (ADU) at location (sourceRecord.get("x"), sourceRecord.get("y")) (pixels) to AB magnitude. More...
 
ndarray::Array< double, 2, 2 > instFluxToMagnitude (afw::table::SourceCatalog const &sourceCatalog, std::string const &instFluxField) const
 Convert sourceCatalog[instFluxField_instFlux] (ADU) at locations (sourceCatalog.get("x"), sourceCatalog.get("y")) (pixels) to AB magnitudes. More...
 
void instFluxToMagnitude (afw::table::SourceCatalog &sourceCatalog, std::string const &instFluxField, std::string const &outField) const
 Convert instFluxes in a catalog to AB magnitudes and write back into the catalog. More...
 
MaskedImage< float > calibrateImage (MaskedImage< float > const &maskedImage, bool includeScaleUncertainty=true) const
 Return a flux calibrated image, with pixel values in nJy. More...
 
afw::table::SourceCatalog calibrateCatalog (afw::table::SourceCatalog const &catalog, std::vector< std::string > const &instFluxFields) const
 Return a flux calibrated catalog, with new _flux, _fluxErr, _mag, and _magErr fields. More...
 
afw::table::SourceCatalog calibrateCatalog (afw::table::SourceCatalog const &catalog) const
 
double magnitudeToInstFlux (double magnitude, lsst::geom::Point< double, 2 > const &point) const
 Convert AB magnitude to instFlux (ADU). More...
 
double magnitudeToInstFlux (double magnitude) const
 
double getCalibrationMean () const
 Get the mean photometric calibration. More...
 
double getCalibrationErr () const
 Get the mean photometric calibration error. More...
 
double getInstFluxAtZeroMagnitude () const
 Get the magnitude zero point (the instrumental flux corresponding to 0 magnitude). More...
 
double getLocalCalibration (lsst::geom::Point< double, 2 > const &point) const
 Get the local calibration at a point. More...
 
std::shared_ptr< afw::math::BoundedFieldcomputeScaledCalibration () const
 Calculates the spatially-variable calibration, normalized by the mean in the valid domain. More...
 
std::shared_ptr< afw::math::BoundedFieldcomputeScalingTo (std::shared_ptr< PhotoCalib > other) const
 Calculates the scaling between this PhotoCalib and another PhotoCalib. More...
 
bool operator== (PhotoCalib const &rhs) const
 Two PhotoCalibs are equal if their component bounded fields and calibrationErr are equal. More...
 
bool operator!= (PhotoCalib const &rhs) const
 Two PhotoCalibs are equal if their component bounded fields and calibrationErr are equal. More...
 
bool isPersistable () const noexcept override
 Return true if this particular object can be persisted using afw::table::io. More...
 
std::shared_ptr< typehandling::StorablecloneStorable () const override
 Create a new PhotoCalib that is a copy of this one. More...
 
std::string toString () const override
 Create a string representation of this object. More...
 
bool equals (typehandling::Storable const &other) const noexcept override
 Compare this object to another Storable. More...
 
virtual std::size_t hash_value () const
 Return a hash of this object (optional operation). More...
 
void writeFits (std::string const &fileName, std::string const &mode="w") const
 Write the object to a regular FITS file. More...
 
void writeFits (fits::MemFileManager &manager, std::string const &mode="w") const
 Write the object to a FITS image in memory. More...
 
void writeFits (fits::Fits &fitsfile) const
 Write the object to an already-open FITS object. More...
 

Static Public Member Functions

static std::shared_ptr< PhotoCalibreadFits (fits::Fits &fitsfile)
 Read an object from an already open FITS object. More...
 
static std::shared_ptr< PhotoCalibreadFits (std::string const &fileName, int hdu=fits::DEFAULT_HDU)
 Read an object from a regular FITS file. More...
 
static std::shared_ptr< PhotoCalibreadFits (fits::MemFileManager &manager, int hdu=fits::DEFAULT_HDU)
 Read an object from a FITS file in memory. More...
 
static std::shared_ptr< PhotoCalibdynamicCast (std::shared_ptr< Persistable > const &ptr)
 Dynamically cast a shared_ptr. More...
 

Protected Types

using OutputArchiveHandle = io::OutputArchiveHandle
 

Protected Member Functions

std::string getPersistenceName () const override
 Return the unique name used to persist this object and look up its factory. More...
 
void write (OutputArchiveHandle &handle) const override
 Write the object to one or more catalogs. More...
 
virtual std::string getPythonModule () const
 Return the fully-qualified Python module that should be imported to guarantee that its factory is registered. More...
 

Static Protected Member Functions

template<class T >
static bool singleClassEquals (T const &lhs, Storable const &rhs)
 Test if a Storable is of a particular class and equal to another object. More...
 

Related Functions

(Note that these are not member functions.)

std::ostreamoperator<< (std::ostream &os, Storable const &storable)
 Output operator for Storable. More...
 

Detailed Description

The photometric calibration of an exposure.

A PhotoCalib is a BoundedField (a function with a specified domain) that converts from post-ISR counts-on-chip (ADU) to flux and magnitude. It is defined such that a calibration of 1 means one count is equal to one nanojansky (nJy, 10^-35 W/m^2/Hz in SI units). The nJy was chosen because it represents a linear flux unit with values in a convenient range (e.g. LSST's single image depth of 24.5 is 575 nJy). See more detailed discussion in: https://pstn-001.lsst.io/

PhotoCalib is immutable.

The spatially varying flux calibration has units of nJy/ADU, and is defined such that, at a position (x,y) in the domain of the boundedField calibration and for a given measured source instFlux:

\[ instFlux*calibration(x,y) = flux [nJy] \]

while the errors (constant on the domain) are defined as:

\[ sqrt((instFluxErr/instFlux)^2 + (calibrationErr/calibration)^2)*flux = fluxErr [nJy] \]

This implies that the conversions from instFlux and instFlux error to magnitude and magnitude error are as follows:

\[ -2.5*log_{10}(instFlux*calibration(x,y)*1e-9/referenceFlux) = magnitude \]

where referenceFlux is the AB Magnitude reference flux from Oke & Gunn 1983 (first equation),

\[ referenceFlux = 1e23 * 10^{(48.6/-2.5)} \]

and

\[ 2.5/log(10)*sqrt((instFluxErr/instFlux)^2 + (calibrationErr/calibration)^2) = magnitudeErr \]

Note that this is independent of referenceFlux.

Definition at line 114 of file PhotoCalib.h.

Member Typedef Documentation

◆ OutputArchiveHandle

using lsst::afw::table::io::Persistable::OutputArchiveHandle = io::OutputArchiveHandle
protectedinherited

Definition at line 108 of file Persistable.h.

Constructor & Destructor Documentation

◆ PhotoCalib() [1/6]

lsst::afw::image::PhotoCalib::PhotoCalib ( PhotoCalib const &  )
default

◆ PhotoCalib() [2/6]

lsst::afw::image::PhotoCalib::PhotoCalib ( PhotoCalib &&  )
default

◆ ~PhotoCalib()

lsst::afw::image::PhotoCalib::~PhotoCalib ( )
overridedefault

◆ PhotoCalib() [3/6]

lsst::afw::image::PhotoCalib::PhotoCalib ( )
inline

Create a empty, zeroed calibration.

Definition at line 127 of file PhotoCalib.h.

127 : PhotoCalib(0) {}
PhotoCalib()
Create a empty, zeroed calibration.
Definition: PhotoCalib.h:127

◆ PhotoCalib() [4/6]

lsst::afw::image::PhotoCalib::PhotoCalib ( double  calibrationMean,
double  calibrationErr = 0,
lsst::geom::Box2I const &  bbox = lsst::geom::Box2I() 
)
inlineexplicit

Create a non-spatially-varying calibration.

Parameters
[in]calibrationMeanThe spatially-constant calibration (must be non-negative).
[in]calibrationErrThe error on the calibration (must be non-negative).
[in]bboxThe bounding box on which this PhotoCalib is valid. If not specified, this PhotoCalib is valid at any point (i.e. an empty bbox).

Definition at line 137 of file PhotoCalib.h.

139  : _calibrationMean(calibrationMean), _calibrationErr(calibrationErr), _isConstant(true) {
140  assertNonNegative(_calibrationMean, "Calibration mean");
141  assertNonNegative(_calibrationErr, "Calibration error");
142  ndarray::Array<double, 2, 2> coeffs = ndarray::allocate(ndarray::makeVector(1, 1));
143  coeffs[0][0] = calibrationMean;
144  _calibration = std::make_shared<afw::math::ChebyshevBoundedField>(
145  afw::math::ChebyshevBoundedField(bbox, coeffs));
146  }
AmpInfoBoxKey bbox
Definition: Amplifier.cc:117
table::Key< double > calibrationMean
Definition: PhotoCalib.cc:369
table::Key< double > calibrationErr
Definition: PhotoCalib.cc:370
void assertNonNegative(double value, std::string const &name)
Raise lsst::pex::exceptions::InvalidParameterError if value is not >=0.
Definition: PhotoCalib.h:69

◆ PhotoCalib() [5/6]

lsst::afw::image::PhotoCalib::PhotoCalib ( std::shared_ptr< afw::math::BoundedField calibration,
double  calibrationErr = 0 
)
inline

Create a spatially-varying calibration.

Parameters
[in]calibrationThe spatially varying photometric calibration (must have non-negative mean).
[in]calibrationErrThe error on the calibration (must be non-negative).

Definition at line 154 of file PhotoCalib.h.

155  : _calibration(calibration),
156  _calibrationMean(computeCalibrationMean(calibration)),
157  _calibrationErr(calibrationErr),
158  _isConstant(false) {
159  assertNonNegative(_calibrationMean, "Calibration (computed via BoundedField.mean()) mean");
160  assertNonNegative(_calibrationErr, "Calibration error");
161  }

◆ PhotoCalib() [6/6]

lsst::afw::image::PhotoCalib::PhotoCalib ( double  calibrationMean,
double  calibrationErr,
std::shared_ptr< afw::math::BoundedField calibration,
bool  isConstant 
)
inline

Create a calibration with a pre-computed mean.

Primarily for de-persistence.

Parameters
[in]calibrationMeanThe mean of the calibration() over its bounding box (must be non-negative).
[in]calibrationErrThe error on the calibration (must be non-negative).
[in]calibrationThe spatially varying photometric calibration.
[in]isConstantIs this PhotoCalib spatially constant?

Definition at line 171 of file PhotoCalib.h.

173  : _calibration(calibration),
174  _calibrationMean(calibrationMean),
175  _calibrationErr(calibrationErr),
176  _isConstant(isConstant) {
177  assertNonNegative(_calibrationMean, "Calibration mean");
178  assertNonNegative(_calibrationErr, "Calibration error");
179  }
table::Key< table::Flag > isConstant
Definition: PhotoCalib.cc:371

Member Function Documentation

◆ calibrateCatalog() [1/2]

afw::table::SourceCatalog lsst::afw::image::PhotoCalib::calibrateCatalog ( afw::table::SourceCatalog const &  catalog) const

Definition at line 349 of file PhotoCalib.cc.

349  {
350  std::vector<std::string> instFluxFields;
351  static std::string const SUFFIX = "_instFlux";
352  for (auto const &name : catalog.getSchema().getNames()) {
353  // Pick every field ending in "_instFlux", grabbing everything before that prefix.
354  if (name.size() > SUFFIX.size() + 1 &&
355  name.compare(name.size() - SUFFIX.size(), SUFFIX.size(), SUFFIX) == 0) {
356  instFluxFields.emplace_back(name.substr(0, name.size() - 9));
357  }
358  }
359  return calibrateCatalog(catalog, instFluxFields);
360 }
table::Key< std::string > name
Definition: Amplifier.cc:116
afw::table::SourceCatalog calibrateCatalog(afw::table::SourceCatalog const &catalog, std::vector< std::string > const &instFluxFields) const
Return a flux calibrated catalog, with new _flux, _fluxErr, _mag, and _magErr fields.
Definition: PhotoCalib.cc:282
T emplace_back(T... args)
T size(T... args)

◆ calibrateCatalog() [2/2]

afw::table::SourceCatalog lsst::afw::image::PhotoCalib::calibrateCatalog ( afw::table::SourceCatalog const &  catalog,
std::vector< std::string > const &  instFluxFields 
) const

Return a flux calibrated catalog, with new _flux, _fluxErr, _mag, and _magErr fields.

If the input catalog already has _flux, _mag, _fluxErr, and/or _magErr fields matching instFluxFields, they will be replaced with the new fields.

Parameters
catalogThe source catalog to compute calibrated fluxes for.
instFluxFieldsThe fields to calibrate (optional). If not provided, every field ending with _instFlux will be calibrated.
Returns
A deep copy of the input catalog, with new calibrated flux and magnitude fields. Elements of instFluxFields that have _instFluxErr will be used to compute the _fluxErr and _magErr fields.
Exceptions
lsst::pex::exceptions::NotFoundErrorif any item in instFluxFields does not have a corresponding *_instFlux field in catalog.schema.

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 282 of file PhotoCalib.cc.

283  {
284  auto const &inSchema = catalog.getSchema();
285  afw::table::SchemaMapper mapper(inSchema, true); // true: share the alias map
286  mapper.addMinimalSchema(inSchema);
287 
288  using FieldD = afw::table::Field<double>;
289 
290  struct Keys {
291  table::Key<double> instFlux;
292  table::Key<double> instFluxErr;
293  table::Key<double> flux;
294  table::Key<double> fluxErr;
295  table::Key<double> mag;
296  table::Key<double> magErr;
297  };
298 
300  keys.reserve(instFluxFields.size());
301  for (auto const &field : instFluxFields) {
302  Keys newKey;
303  newKey.instFlux = inSchema[inSchema.join(field, "instFlux")];
304  newKey.flux =
305  mapper.addOutputField(FieldD(inSchema.join(field, "flux"), "calibrated flux", "nJy"), true);
306  newKey.mag = mapper.addOutputField(
307  FieldD(inSchema.join(field, "mag"), "calibrated magnitude", "mag(AB)"), true);
308  try {
309  newKey.instFluxErr = inSchema.find<double>(inSchema.join(field, "instFluxErr")).key;
310  newKey.fluxErr = mapper.addOutputField(
311  FieldD(inSchema.join(field, "fluxErr"), "calibrated flux uncertainty", "nJy"), true);
312  newKey.magErr = mapper.addOutputField(
313  FieldD(inSchema.join(field, "magErr"), "calibrated magnitude uncertainty", "mag(AB)"),
314  true);
315  } catch (pex::exceptions::NotFoundError &) {
316  ; // Keys struct defaults to invalid keys; that marks the error as missing.
317  }
318  keys.emplace_back(newKey);
319  }
320 
321  // Create the new catalog
322  afw::table::SourceCatalog output(mapper.getOutputSchema());
323  output.insert(mapper, output.begin(), catalog.begin(), catalog.end());
324 
325  auto calibration = evaluateCatalog(output);
326 
327  // fill in the catalog values
328  int iRec = 0;
329  for (auto &rec : output) {
330  for (auto &key : keys) {
331  double instFlux = rec.get(key.instFlux);
332  double nanojansky = toNanojansky(instFlux, calibration[iRec]);
333  rec.set(key.flux, nanojansky);
334  rec.set(key.mag, toMagnitude(instFlux, calibration[iRec]));
335  if (key.instFluxErr.isValid()) {
336  double instFluxErr = rec.get(key.instFluxErr);
337  rec.set(key.fluxErr, toNanojanskyErr(instFlux, instFluxErr, calibration[iRec],
338  _calibrationErr, nanojansky));
339  rec.set(key.magErr,
340  toMagnitudeErr(instFlux, instFluxErr, calibration[iRec], _calibrationErr));
341  }
342  }
343  ++iRec;
344  }
345 
346  return output;
347 }
table::Key< int > field
Definition: PhotoCalib.cc:372
SchemaMapper * mapper
Definition: SchemaMapper.cc:71
SortedCatalogT< SourceRecord > SourceCatalog
Definition: fwd.h:85

◆ calibrateImage()

MaskedImage< float > lsst::afw::image::PhotoCalib::calibrateImage ( MaskedImage< float > const &  maskedImage,
bool  includeScaleUncertainty = true 
) const

Return a flux calibrated image, with pixel values in nJy.

Mask pixels are propagated directly from the input image.

Parameters
maskedImageThe masked image to calibrate.
includeScaleUncertaintyInclude the uncertainty on the calibration in the resulting variance?
Returns
The calibrated masked image.

Definition at line 260 of file PhotoCalib.cc.

261  {
262  // Deep copy construct, as we're mutiplying in-place.
263  auto result = MaskedImage<float>(maskedImage, true);
264 
265  if (_isConstant) {
266  *(result.getImage()) *= _calibrationMean;
267  } else {
268  _calibration->multiplyImage(*(result.getImage()), true); // only in the overlap region
269  }
270  if (includeScaleUncertainty) {
271  toNanojanskyVariance(maskedImage.getImage()->getArray(), maskedImage.getVariance()->getArray(),
272  _calibrationErr, result.getImage()->getArray(),
273  result.getVariance()->getArray());
274  } else {
275  toNanojanskyVariance(maskedImage.getImage()->getArray(), maskedImage.getVariance()->getArray(), 0,
276  result.getImage()->getArray(), result.getVariance()->getArray());
277  }
278 
279  return result;
280 }
py::object result
Definition: _schema.cc:429

◆ cloneStorable()

std::shared_ptr< typehandling::Storable > lsst::afw::image::PhotoCalib::cloneStorable ( ) const
overridevirtual

Create a new PhotoCalib that is a copy of this one.

Reimplemented from lsst::afw::typehandling::Storable.

Definition at line 238 of file PhotoCalib.cc.

238  {
239  return std::make_unique<PhotoCalib>(*this);
240 }

◆ computeScaledCalibration()

std::shared_ptr< math::BoundedField > lsst::afw::image::PhotoCalib::computeScaledCalibration ( ) const

Calculates the spatially-variable calibration, normalized by the mean in the valid domain.

This value is defined, for instFlux at (x,y), such that:

\[ instFlux*computeScaledCalibration()(x,y)*getCalibrationMean() = instFluxToNanojansky(instFlux, (x,y)) \]

See also
PhotoCalib::getCalibrationMean()
Returns
The normalized spatially-variable calibration.

Definition at line 221 of file PhotoCalib.cc.

221  {
222  return *(_calibration) / _calibrationMean;
223 }

◆ computeScalingTo()

std::shared_ptr< math::BoundedField > lsst::afw::image::PhotoCalib::computeScalingTo ( std::shared_ptr< PhotoCalib other) const

Calculates the scaling between this PhotoCalib and another PhotoCalib.

The BoundedFields of these PhotoCalibs must have the same BBoxes (or one or both must be empty).

With:

  • c = instFlux at position (x,y)
  • this = this PhotoCalib
  • other = other PhotoCalib
  • return = BoundedField returned by this method the return value from this method is defined as:

    \[ this.instFluxToNanojansky(c, (x,y))*return(x, y) = other.instFluxToNanojansky(c, (x,y)) \]

Parameters
[in]otherThe PhotoCalib to scale to.
Returns
The BoundedField as defined above.
Warning
Not implemented yet: See DM-10154.

Definition at line 225 of file PhotoCalib.cc.

225  {
226  throw LSST_EXCEPT(pex::exceptions::LogicError, "Not Implemented: See DM-10154.");
227 }
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48

◆ dynamicCast()

Dynamically cast a shared_ptr.

Dynamically cast a shared pointer and raise on failure.

You must provide an explicit template instantiation in the .cc file for each class that inherits from PersistableFacade. Designed to work around RTTI issues on macOS with hidden symbols;

Exceptions
lsst::pex::exceptions::LogicErrorif the cast fails

param[in] ptr The pointer to be cast.

Returns
The cast pointer.
Exceptions
lsst::pex::exceptions::TypeErrorIf the dynamic cast fails.

Definition at line 218 of file Persistable.cc.

18  {
19  auto result = std::dynamic_pointer_cast<T>(ptr);
20  if (!result) {
21  throw LSST_EXCEPT(pex::exceptions::TypeError, "Dynamic pointer cast failed");
22  }
23  return result;
24 }
uint64_t * ptr
Definition: RangeSet.cc:88

◆ equals()

bool lsst::afw::image::PhotoCalib::equals ( typehandling::Storable const &  other) const
overridevirtualnoexcept

Compare this object to another Storable.

Returns
*this == other if other is a PhotoCalib; otherwise false.

Reimplemented from lsst::afw::typehandling::Storable.

Definition at line 252 of file PhotoCalib.cc.

252  {
253  return singleClassEquals(*this, other);
254 }
static bool singleClassEquals(T const &lhs, Storable const &rhs)
Test if a Storable is of a particular class and equal to another object.
Definition: Storable.h:151

◆ getCalibrationErr()

double lsst::afw::image::PhotoCalib::getCalibrationErr ( ) const
inline

Get the mean photometric calibration error.

This value is defined such that for some instFluxErr, instFlux, and flux:

\[ sqrt((instFluxErr/instFlux)^2 + (calibrationErr/calibration(x,y))^2)*flux = fluxErr [nJy] \]

See also
PhotoCalib::computeScaledCalibration(), getCalibrationMean()
Returns
The calibration error.

Definition at line 418 of file PhotoCalib.h.

418 { return _calibrationErr; }

◆ getCalibrationMean()

double lsst::afw::image::PhotoCalib::getCalibrationMean ( ) const
inline

Get the mean photometric calibration.

This value is defined, for instFlux at (x,y), such that:

\[ instFlux*computeScaledCalibration()(x,y)*getCalibrationMean() = instFluxToNanojansky(instFlux, (x,y)) \]

See also
PhotoCalib::computeScaledCalibration(), getCalibrationErr(), getInstFluxAtZeroMagnitude()
Returns
The spatial mean of this calibration.

Definition at line 404 of file PhotoCalib.h.

404 { return _calibrationMean; }

◆ getInstFluxAtZeroMagnitude()

double lsst::afw::image::PhotoCalib::getInstFluxAtZeroMagnitude ( ) const
inline

Get the magnitude zero point (the instrumental flux corresponding to 0 magnitude).

This value is defined such that:

\[ instFluxToMagnitude(getInstFluxAtZeroMagnitude()) == 0 \]

See also
PhotoCalib::computeScaledCalibration(), getCalibrationMean()
Returns
The instFlux magnitude zero point.

Definition at line 432 of file PhotoCalib.h.

432 { return utils::referenceFlux / _calibrationMean; }
const double referenceFlux
The Oke & Gunn (1983) AB magnitude reference flux, in nJy (often approximated as 3631....
Definition: Magnitude.h:46

◆ getLocalCalibration()

double lsst::afw::image::PhotoCalib::getLocalCalibration ( lsst::geom::Point< double, 2 > const &  point) const
inline

Get the local calibration at a point.

This value is defined such that:

\[ instFluxToNanojansky(instFlux, point) == localCalibration * inst \]

Use getCalibrationErr() to get the spatially constant error on the calibration.

Parameters
[in]pointThe position that magnitude is to be converted at.
See also
getCalibrationMean(), getCalibrationErr()
Returns
The local calibration at a point.

Definition at line 450 of file PhotoCalib.h.

450 { return evaluate(point); }

◆ getPersistenceName()

std::string lsst::afw::image::PhotoCalib::getPersistenceName ( ) const
overrideprotectedvirtual

Return the unique name used to persist this object and look up its factory.

Must be less than ArchiveIndexSchema::MAX_NAME_LENGTH characters.

Reimplemented from lsst::afw::table::io::Persistable.

Definition at line 509 of file PhotoCalib.cc.

509 { return getPhotoCalibPersistenceName(); }

◆ getPythonModule()

std::string lsst::afw::table::io::Persistable::getPythonModule ( ) const
protectedvirtualinherited

◆ hash_value()

std::size_t lsst::afw::typehandling::Storable::hash_value ( ) const
virtualinherited

Return a hash of this object (optional operation).

Exceptions
UnsupportedOperationExceptionThrown if this object is not hashable.
Note
C++ subclass authors are responsible for any associated specializations of std::hash.
When called on Python classes, this method delegates to __hash__ if it exists.

Reimplemented in lsst::afw::typehandling::StorableHelper< Base >, lsst::afw::image::VisitInfo, lsst::afw::image::FilterLabel, and lsst::afw::geom::polygon::Polygon.

Definition at line 44 of file Storable.cc.

44  {
45  throw LSST_EXCEPT(UnsupportedOperationException, "Hashes are not supported.");
46 }

◆ instFluxToMagnitude() [1/7]

void lsst::afw::image::PhotoCalib::instFluxToMagnitude ( afw::table::SourceCatalog sourceCatalog,
std::string const &  instFluxField,
std::string const &  outField 
) const

Convert instFluxes in a catalog to AB magnitudes and write back into the catalog.

Convert sourceCatalog[instFluxField_instFlux] (ADU) at locations (sourceCatalog.get("x"), sourceCatalog.get("y")) (pixels) to AB magnitudes and write the results back to sourceCatalog[outField_mag].

Parameters
[in]sourceCatalogThe source catalog to get instFlux and position from.
[in]instFluxFieldThe instFlux field: Keys of the form "*_instFlux" and "*_instFluxErr" must exist. For example: instFluxField="slot_PsfFlux" will use the fields named: "slot_PsfFlux_instFlux", "slot_PsfFlux_instFluxErr"
[in]outFieldThe field to write the magnitudes and magnitude errors to. Keys of the form "*_instFlux", "*_instFluxErr", *_mag", and "*_magErr" must exist in the schema.

Definition at line 197 of file PhotoCalib.cc.

198  {
199  auto instFluxKey = sourceCatalog.getSchema().find<double>(instFluxField + "_instFlux").key;
200  auto instFluxErrKey = sourceCatalog.getSchema().find<double>(instFluxField + "_instFluxErr").key;
201  auto magKey = sourceCatalog.getSchema().find<double>(outField + "_mag").key;
202  auto magErrKey = sourceCatalog.getSchema().find<double>(outField + "_magErr").key;
203  for (auto &record : sourceCatalog) {
204  auto result = instFluxToMagnitude(record.get(instFluxKey), record.get(instFluxErrKey),
205  record.getCentroid());
206  record.set(magKey, result.value);
207  record.set(magErrKey, result.error);
208  }
209 }
double instFluxToMagnitude(double instFlux, lsst::geom::Point< double, 2 > const &point) const
Convert instFlux in ADU to AB magnitude.
Definition: PhotoCalib.cc:158

◆ instFluxToMagnitude() [2/7]

ndarray::Array< double, 2, 2 > lsst::afw::image::PhotoCalib::instFluxToMagnitude ( afw::table::SourceCatalog const &  sourceCatalog,
std::string const &  instFluxField 
) const

Convert sourceCatalog[instFluxField_instFlux] (ADU) at locations (sourceCatalog.get("x"), sourceCatalog.get("y")) (pixels) to AB magnitudes.

Parameters
[in]sourceCatalogThe source catalog to get instFlux and position from.
[in]instFluxFieldThe instFlux field: Keys of the form "*_instFlux" and "*_instFluxErr" must exist. For example: instFluxField="slot_PsfFlux" will use the fields named: "slot_PsfFlux_instFlux", "slot_PsfFlux_instFluxErr"
Returns
The magnitudes and magnitude errors for the sources.

Definition at line 189 of file PhotoCalib.cc.

190  {
191  ndarray::Array<double, 2, 2> result =
192  ndarray::allocate(ndarray::makeVector(int(sourceCatalog.size()), 2));
193  instFluxToMagnitudeArray(sourceCatalog, instFluxField, result);
194  return result;
195 }

◆ instFluxToMagnitude() [3/7]

Measurement lsst::afw::image::PhotoCalib::instFluxToMagnitude ( afw::table::SourceRecord const &  sourceRecord,
std::string const &  instFluxField 
) const

Convert sourceRecord[instFluxField_instFlux] (ADU) at location (sourceRecord.get("x"), sourceRecord.get("y")) (pixels) to AB magnitude.

Parameters
[in]sourceRecordThe source record to get instFlux and position from.
[in]instFluxFieldThe instFlux field: Keys of the form "*_instFlux" and "*_instFluxErr" must exist. For example: instFluxField="slot_PsfFlux" will use the fields named: "slot_PsfFlux_instFlux", "slot_PsfFlux_instFluxErr"
Returns
The magnitude and magnitude error for this source.

Definition at line 181 of file PhotoCalib.cc.

182  {
183  auto position = sourceRecord.getCentroid();
184  auto instFluxKey = sourceRecord.getSchema().find<double>(instFluxField + "_instFlux").key;
185  auto instFluxErrKey = sourceRecord.getSchema().find<double>(instFluxField + "_instFluxErr").key;
186  return instFluxToMagnitude(sourceRecord.get(instFluxKey), sourceRecord.get(instFluxErrKey), position);
187 }

◆ instFluxToMagnitude() [4/7]

double lsst::afw::image::PhotoCalib::instFluxToMagnitude ( double  instFlux) const

Definition at line 162 of file PhotoCalib.cc.

162  {
163  return toMagnitude(instFlux, _calibrationMean);
164 }

◆ instFluxToMagnitude() [5/7]

Measurement lsst::afw::image::PhotoCalib::instFluxToMagnitude ( double  instFlux,
double  instFluxErr 
) const

Definition at line 175 of file PhotoCalib.cc.

175  {
176  double magnitude = toMagnitude(instFlux, _calibrationMean);
177  double error = toMagnitudeErr(instFlux, instFluxErr, _calibrationMean, _calibrationErr);
178  return Measurement(magnitude, error);
179 }

◆ instFluxToMagnitude() [6/7]

Measurement lsst::afw::image::PhotoCalib::instFluxToMagnitude ( double  instFlux,
double  instFluxErr,
lsst::geom::Point< double, 2 > const &  point 
) const

Convert instFlux and error in instFlux (ADU) to AB magnitude and magnitude error.

If passed point, use the exact calculation at that point, otherwise, use the mean scaling factor.

Parameters
[in]instFluxThe source instFlux in ADU.
[in]instFluxErrThe instFlux error (standard deviation).
[in]pointThe point that instFlux is measured at.
Returns
The AB magnitude and error.

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 166 of file PhotoCalib.cc.

167  {
168  double calibration, error, magnitude;
169  calibration = evaluate(point);
170  magnitude = toMagnitude(instFlux, calibration);
171  error = toMagnitudeErr(instFlux, instFluxErr, calibration, _calibrationErr);
172  return Measurement(magnitude, error);
173 }

◆ instFluxToMagnitude() [7/7]

double lsst::afw::image::PhotoCalib::instFluxToMagnitude ( double  instFlux,
lsst::geom::Point< double, 2 > const &  point 
) const

Convert instFlux in ADU to AB magnitude.

If passed point, use the exact calculation at that point, otherwise, use the mean scaling factor.

Parameters
[in]instFluxThe source instFlux in ADU.
[in]pointThe point that instFlux is measured at.
Returns
The AB magnitude.

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 158 of file PhotoCalib.cc.

158  {
159  return toMagnitude(instFlux, evaluate(point));
160 }

◆ instFluxToNanojansky() [1/7]

void lsst::afw::image::PhotoCalib::instFluxToNanojansky ( afw::table::SourceCatalog sourceCatalog,
std::string const &  instFluxField,
std::string const &  outField 
) const

Convert sourceCatalog[instFluxField_instFlux] (ADU) at locations (sourceCatalog.get("x"), sourceCatalog.get("y")) (pixels) to nJy and write the results back to sourceCatalog[outField_mag].

Parameters
[in]sourceCatalogThe source catalog to get instFlux and position from.
[in]instFluxFieldThe instFlux field: Keys of the form "*_instFlux" and "*_instFluxErr" must exist. For example: instFluxField="slot_PsfFlux" will use the fields named: "slot_PsfFlux_instFlux", "slot_PsfFlux_instFluxErr"
[in]outFieldThe field to write the nJy and magnitude errors to. Keys of the form "*_instFlux" and "*_instFluxErr" must exist in the schema.

Definition at line 142 of file PhotoCalib.cc.

143  {
144  auto instFluxKey = sourceCatalog.getSchema().find<double>(instFluxField + "_instFlux").key;
145  auto instFluxErrKey = sourceCatalog.getSchema().find<double>(instFluxField + "_instFluxErr").key;
146  auto nanojanskyKey = sourceCatalog.getSchema().find<double>(outField + "_flux").key;
147  auto nanojanskyErrKey = sourceCatalog.getSchema().find<double>(outField + "_fluxErr").key;
148  for (auto &record : sourceCatalog) {
149  auto result = instFluxToNanojansky(record.get(instFluxKey), record.get(instFluxErrKey),
150  record.getCentroid());
151  record.set(nanojanskyKey, result.value);
152  record.set(nanojanskyErrKey, result.error);
153  }
154 }
double instFluxToNanojansky(double instFlux, lsst::geom::Point< double, 2 > const &point) const
Convert instFlux in ADU to nJy at a point in the BoundedField.
Definition: PhotoCalib.cc:104

◆ instFluxToNanojansky() [2/7]

ndarray::Array< double, 2, 2 > lsst::afw::image::PhotoCalib::instFluxToNanojansky ( afw::table::SourceCatalog const &  sourceCatalog,
std::string const &  instFluxField 
) const

Convert sourceCatalog[instFluxField_instFlux] (ADU) at locations (sourceCatalog.get("x"), sourceCatalog.get("y")) (pixels) to nJy.

Parameters
[in]sourceCatalogThe source catalog to get instFlux and position from.
[in]instFluxFieldThe instFlux field: Keys of the form "*_instFlux" and "*_instFluxErr" must exist. For example: instFluxField="slot_PsfFlux" will use the fields named: "slot_PsfFlux_instFlux", "slot_PsfFlux_instFluxErr"
Returns
The flux in nJy and error for this source.

Definition at line 134 of file PhotoCalib.cc.

135  {
136  ndarray::Array<double, 2, 2> result =
137  ndarray::allocate(ndarray::makeVector(int(sourceCatalog.size()), 2));
138  instFluxToNanojanskyArray(sourceCatalog, instFluxField, result);
139  return result;
140 }

◆ instFluxToNanojansky() [3/7]

Measurement lsst::afw::image::PhotoCalib::instFluxToNanojansky ( const afw::table::SourceRecord sourceRecord,
std::string const &  instFluxField 
) const

Convert sourceRecord[instFluxField_instFlux] (ADU) at location (sourceRecord.get("x"), sourceRecord.get("y")) (pixels) to flux and flux error (in nJy).

Parameters
[in]sourceRecordThe source record to get instFlux and position from.
[in]instFluxFieldThe instFlux field: Keys of the form "*_instFlux" and "*_instFluxErr" must exist. For example: instFluxField="slot_PsfFlux" will use the fields named: "slot_PsfFlux_instFlux", "slot_PsfFlux_instFluxErr"
Returns
The flux in nJy and error for this source.

Definition at line 127 of file PhotoCalib.cc.

128  {
129  auto position = sourceRecord.getCentroid();
130  auto instFluxKey = sourceRecord.getSchema().find<double>(instFluxField + "_instFlux").key;
131  auto instFluxErrKey = sourceRecord.getSchema().find<double>(instFluxField + "_instFluxErr").key;
132  return instFluxToNanojansky(sourceRecord.get(instFluxKey), sourceRecord.get(instFluxErrKey), position);
133 }

◆ instFluxToNanojansky() [4/7]

double lsst::afw::image::PhotoCalib::instFluxToNanojansky ( double  instFlux) const

Definition at line 108 of file PhotoCalib.cc.

108  {
109  return toNanojansky(instFlux, _calibrationMean);
110 }

◆ instFluxToNanojansky() [5/7]

Measurement lsst::afw::image::PhotoCalib::instFluxToNanojansky ( double  instFlux,
double  instFluxErr 
) const

Definition at line 121 of file PhotoCalib.cc.

121  {
122  double nanojansky = toNanojansky(instFlux, _calibrationMean);
123  double error = toNanojanskyErr(instFlux, instFluxErr, _calibrationMean, _calibrationErr, nanojansky);
124  return Measurement(nanojansky, error);
125 }

◆ instFluxToNanojansky() [6/7]

Measurement lsst::afw::image::PhotoCalib::instFluxToNanojansky ( double  instFlux,
double  instFluxErr,
lsst::geom::Point< double, 2 > const &  point 
) const

Convert instFlux and error in instFlux (ADU) to nJy and nJy error.

If passed point, use the exact calculation at that point, otherwise, use the mean scaling factor.

Parameters
[in]instFluxThe source fluxinstFlux in ADU.
[in]instFluxErrThe instFlux error.
[in]pointThe point that instFlux is measured at.
Returns
The flux in nJy and error.

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 112 of file PhotoCalib.cc.

113  {
114  double calibration, error, nanojansky;
115  calibration = evaluate(point);
116  nanojansky = toNanojansky(instFlux, calibration);
117  error = toNanojanskyErr(instFlux, instFluxErr, calibration, _calibrationErr, nanojansky);
118  return Measurement(nanojansky, error);
119 }

◆ instFluxToNanojansky() [7/7]

double lsst::afw::image::PhotoCalib::instFluxToNanojansky ( double  instFlux,
lsst::geom::Point< double, 2 > const &  point 
) const

Convert instFlux in ADU to nJy at a point in the BoundedField.

If passed point, use the exact calculation at that point, otherwise, use the mean scaling factor.

Parameters
[in]instFluxThe source instFlux in ADU.
[in]pointThe point that instFlux is measured at.
Returns
The flux in nJy.

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 104 of file PhotoCalib.cc.

104  {
105  return toNanojansky(instFlux, evaluate(point));
106 }

◆ isPersistable()

bool lsst::afw::image::PhotoCalib::isPersistable ( ) const
inlineoverridevirtualnoexcept

Return true if this particular object can be persisted using afw::table::io.

Reimplemented from lsst::afw::table::io::Persistable.

Definition at line 496 of file PhotoCalib.h.

496 { return true; }

◆ magnitudeToInstFlux() [1/2]

double lsst::afw::image::PhotoCalib::magnitudeToInstFlux ( double  magnitude) const

Definition at line 213 of file PhotoCalib.cc.

213  {
214  return toInstFluxFromMagnitude(magnitude, _calibrationMean);
215 }

◆ magnitudeToInstFlux() [2/2]

double lsst::afw::image::PhotoCalib::magnitudeToInstFlux ( double  magnitude,
lsst::geom::Point< double, 2 > const &  point 
) const

Convert AB magnitude to instFlux (ADU).

If passed point, use the exact calculation at that point, otherwise, use the mean scaling factor.

Useful for inserting fake sources into an image.

Parameters
[in]magnitudeThe AB magnitude to convert.
[in]pointThe position that magnitude is to be converted at.
Returns
Source instFlux in ADU.

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Definition at line 217 of file PhotoCalib.cc.

217  {
218  return toInstFluxFromMagnitude(magnitude, evaluate(point));
219 }

◆ operator!=()

bool lsst::afw::image::PhotoCalib::operator!= ( PhotoCalib const &  rhs) const
inline

Two PhotoCalibs are equal if their component bounded fields and calibrationErr are equal.

Definition at line 494 of file PhotoCalib.h.

494 { return !(*this == rhs); }

◆ operator=() [1/2]

PhotoCalib& lsst::afw::image::PhotoCalib::operator= ( PhotoCalib &&  )
delete

◆ operator=() [2/2]

PhotoCalib& lsst::afw::image::PhotoCalib::operator= ( PhotoCalib const &  )
delete

◆ operator==()

bool lsst::afw::image::PhotoCalib::operator== ( PhotoCalib const &  rhs) const

Two PhotoCalibs are equal if their component bounded fields and calibrationErr are equal.

Definition at line 229 of file PhotoCalib.cc.

229  {
230  return (_calibrationMean == rhs._calibrationMean && _calibrationErr == rhs._calibrationErr &&
231  (*_calibration) == *(rhs._calibration));
232 }

◆ readFits() [1/3]

static std::shared_ptr<PhotoCalib > lsst::afw::table::io::PersistableFacade< PhotoCalib >::readFits ( fits::Fits fitsfile)
inlinestaticinherited

Read an object from an already open FITS object.

Parameters
[in]fitsfileFITS object to read from, already positioned at the desired HDU.

Definition at line 183 of file Persistable.h.

183  {
184  return dynamicCast(Persistable::_readFits(fitsfile));
185  }
static std::shared_ptr< PhotoCalib > dynamicCast(std::shared_ptr< Persistable > const &ptr)
Dynamically cast a shared_ptr.
Definition: Persistable.cc:18

◆ readFits() [2/3]

static std::shared_ptr<PhotoCalib > lsst::afw::table::io::PersistableFacade< PhotoCalib >::readFits ( fits::MemFileManager manager,
int  hdu = fits::DEFAULT_HDU 
)
inlinestaticinherited

Read an object from a FITS file in memory.

Parameters
[in]managerManager for the memory to read from.
[in]hduHDU to read, where 0 is the primary. The special value of afw::fits::DEFAULT_HDU skips the primary HDU if it is empty.

Definition at line 205 of file Persistable.h.

205  {
206  return dynamicCast(Persistable::_readFits(manager, hdu));
207  }

◆ readFits() [3/3]

static std::shared_ptr<PhotoCalib > lsst::afw::table::io::PersistableFacade< PhotoCalib >::readFits ( std::string const &  fileName,
int  hdu = fits::DEFAULT_HDU 
)
inlinestaticinherited

Read an object from a regular FITS file.

Parameters
[in]fileNameName of the file to read.
[in]hduHDU to read, where 0 is the primary. The special value of afw::fits::DEFAULT_HDU skips the primary HDU if it is empty.

Definition at line 194 of file Persistable.h.

194  {
195  return dynamicCast(Persistable::_readFits(fileName, hdu));
196  }

◆ singleClassEquals()

template<class T >
static bool lsst::afw::typehandling::Storable::singleClassEquals ( T const &  lhs,
Storable const &  rhs 
)
inlinestaticprotectedinherited

Test if a Storable is of a particular class and equal to another object.

This method template simplifies implementations of equals that delegate to operator== without supporting cross-class comparisons.

Template Parameters
TThe class expected of the two objects to be compared.
Parameters
lhs,rhsThe objects to compare. Note that rhs need not be a T, while lhs must be.
Returns
true if rhs is a T and lhs == rhs; false otherwise.
Exception Safety
Provides the same level of exception safety as operator==. Most implementations of operator== do not throw.
Note
This method template calls operator== with both arguments of compile-time type T const&. Its use is not recommended if there would be any ambiguity as to which operator== gets picked by overload resolution.

This method template is typically called from equals as:

bool MyType::equals(Storable const& other) const noexcept {
    return singleClassEquals(*this, other);
}

Definition at line 151 of file Storable.h.

151  {
152  auto typedRhs = dynamic_cast<T const*>(&rhs);
153  if (typedRhs != nullptr) {
154  return lhs == *typedRhs;
155  } else {
156  return false;
157  }
158  }

◆ toString()

std::string lsst::afw::image::PhotoCalib::toString ( ) const
overridevirtual

Create a string representation of this object.

Reimplemented from lsst::afw::typehandling::Storable.

Definition at line 242 of file PhotoCalib.cc.

242  {
243  std::stringstream buffer;
244  if (_isConstant)
245  buffer << "spatially constant with ";
246  else
247  buffer << *_calibration << " with ";
248  buffer << "mean: " << _calibrationMean << " error: " << _calibrationErr;
249  return buffer.str();
250 }
T str(T... args)

◆ write()

void lsst::afw::image::PhotoCalib::write ( OutputArchiveHandle handle) const
overrideprotectedvirtual

Write the object to one or more catalogs.

The handle object passed to this function provides an interface for adding new catalogs and adding nested objects to the same archive (while checking for duplicates). See OutputArchiveHandle for more information.

Reimplemented from lsst::afw::table::io::Persistable.

Definition at line 511 of file PhotoCalib.cc.

511  {
512  PhotoCalibSchema const &keys = PhotoCalibSchema::get();
513  table::BaseCatalog catalog = handle.makeCatalog(keys.schema);
514  auto record = catalog.addNew();
515  record->set(keys.calibrationMean, _calibrationMean);
516  record->set(keys.calibrationErr, _calibrationErr);
517  record->set(keys.isConstant, _isConstant);
518  record->set(keys.field, handle.put(_calibration));
519  record->set(keys.version, SERIALIZATION_VERSION);
520  handle.saveCatalog(catalog);
521 }
std::shared_ptr< RecordT > addNew()
Create a new record, add it to the end of the catalog, and return a pointer to it.
Definition: Catalog.h:490
CatalogT< BaseRecord > BaseCatalog
Definition: fwd.h:72

◆ writeFits() [1/3]

void lsst::afw::table::io::Persistable::writeFits ( fits::Fits fitsfile) const
inherited

Write the object to an already-open FITS object.

Parameters
[in]fitsfileOpen FITS object to write to.

Definition at line 18 of file Persistable.cc.

18  {
19  OutputArchive archive;
20  archive.put(this);
21  archive.writeFits(fitsfile);
22 }

◆ writeFits() [2/3]

void lsst::afw::table::io::Persistable::writeFits ( fits::MemFileManager manager,
std::string const &  mode = "w" 
) const
inherited

Write the object to a FITS image in memory.

Parameters
[in]managerName of the file to write to.
[in]modeIf "w", any existing file with the given name will be overwritten. If "a", new HDUs will be appended to an existing file.

Definition at line 29 of file Persistable.cc.

29  {
30  fits::Fits fitsfile(manager, mode, fits::Fits::AUTO_CLOSE | fits::Fits::AUTO_CHECK);
31  writeFits(fitsfile);
32 }
void writeFits(std::string const &fileName, std::string const &mode="w") const
Write the object to a regular FITS file.
Definition: Persistable.cc:24

◆ writeFits() [3/3]

void lsst::afw::table::io::Persistable::writeFits ( std::string const &  fileName,
std::string const &  mode = "w" 
) const
inherited

Write the object to a regular FITS file.

Parameters
[in]fileNameName of the file to write to.
[in]modeIf "w", any existing file with the given name will be overwritten. If "a", new HDUs will be appended to an existing file.

Definition at line 24 of file Persistable.cc.

24  {
25  fits::Fits fitsfile(fileName, mode, fits::Fits::AUTO_CLOSE | fits::Fits::AUTO_CHECK);
26  writeFits(fitsfile);
27 }

Friends And Related Function Documentation

◆ operator<<()

std::ostream & operator<< ( std::ostream os,
Storable const &  storable 
)
related

Output operator for Storable.

Parameters
osthe desired output stream
storablethe object to print
Returns
a reference to os
Exceptions
UnsupportedOperationExceptionThrown if storable does not have an implementation of Storable::toString.

Definition at line 174 of file Storable.h.

174  {
175  return os << storable.toString();
176 }
std::ostream * os
Definition: Schema.cc:557

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