LSST Applications g0265f82a02+0e5473021a,g02d81e74bb+bd2ed33bd6,g1470d8bcf6+c6d6eb38e2,g14a832a312+9d12ad093c,g2079a07aa2+86d27d4dc4,g2305ad1205+91a32aca49,g295015adf3+88246b6574,g2bbee38e9b+0e5473021a,g337abbeb29+0e5473021a,g3ddfee87b4+c34e8be1fa,g487adcacf7+85dcfbcc36,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+ea1711114f,g5a732f18d5+53520f316c,g64a986408d+bd2ed33bd6,g858d7b2824+bd2ed33bd6,g8a8a8dda67+585e252eca,g99cad8db69+016a06b37a,g9ddcbc5298+9a081db1e4,ga1e77700b3+15fc3df1f7,ga8c6da7877+ef4e3a5875,gb0e22166c9+60f28cb32d,gb6a65358fc+0e5473021a,gba4ed39666+c2a2e4ac27,gbb8dafda3b+b6d7b42999,gc120e1dc64+f745648b3a,gc28159a63d+0e5473021a,gcf0d15dbbd+c34e8be1fa,gdaeeff99f8+f9a426f77a,ge6526c86ff+508d0e0a30,ge79ae78c31+0e5473021a,gee10cc3b42+585e252eca,gf18bd8381d+8d59551888,gf1cff7945b+bd2ed33bd6,w.2024.16
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | Static Protected Member Functions | Related Symbols | 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.
 
 PhotoCalib (double calibrationMean, double calibrationErr=0, lsst::geom::Box2I const &bbox=lsst::geom::Box2I())
 Create a non-spatially-varying calibration.
 
 PhotoCalib (std::shared_ptr< afw::math::BoundedField > calibration, double calibrationErr=0)
 Create a spatially-varying calibration.
 
 PhotoCalib (double calibrationMean, double calibrationErr, std::shared_ptr< afw::math::BoundedField > calibration, bool isConstant)
 Create a calibration with a pre-computed mean.
 
double instFluxToNanojansky (double instFlux, lsst::geom::Point< double, 2 > const &point) const
 Convert instFlux in ADU to nJy at a point in the BoundedField.
 
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.
 
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).
 
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.
 
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].
 
double instFluxToMagnitude (double instFlux, lsst::geom::Point< double, 2 > const &point) const
 Convert instFlux in ADU to AB magnitude.
 
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.
 
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.
 
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.
 
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.
 
MaskedImage< float > calibrateImage (MaskedImage< float > const &maskedImage, bool includeScaleUncertainty=true) const
 Return a flux calibrated image, with pixel values in nJy.
 
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.
 
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).
 
double magnitudeToInstFlux (double magnitude) const
 
double getCalibrationMean () const
 Get the mean photometric calibration.
 
double getCalibrationErr () const
 Get the mean photometric calibration error.
 
bool isConstant () const
 Is this photoCalib spatially constant?
 
double getInstFluxAtZeroMagnitude () const
 Get the magnitude zero point (the instrumental flux corresponding to 0 magnitude).
 
double getLocalCalibration (lsst::geom::Point< double, 2 > const &point) const
 Get the local calibration at a point.
 
std::shared_ptr< afw::math::BoundedFieldcomputeScaledCalibration () const
 Calculates the spatially-variable calibration, normalized by the mean in the valid domain.
 
std::shared_ptr< afw::math::BoundedFieldcomputeScalingTo (std::shared_ptr< PhotoCalib > other) const
 Calculates the scaling between this PhotoCalib and another PhotoCalib.
 
bool operator== (PhotoCalib const &rhs) const
 Two PhotoCalibs are equal if their component bounded fields and calibrationErr are equal.
 
bool operator!= (PhotoCalib const &rhs) const
 Two PhotoCalibs are equal if their component bounded fields and calibrationErr are equal.
 
bool isPersistable () const noexcept override
 Return true if this particular object can be persisted using afw::table::io.
 
std::shared_ptr< typehandling::StorablecloneStorable () const override
 Create a new PhotoCalib that is a copy of this one.
 
std::string toString () const override
 Create a string representation of this object.
 
bool equals (typehandling::Storable const &other) const noexcept override
 Compare this object to another Storable.
 
virtual std::size_t hash_value () const
 Return a hash of this object (optional operation).
 
void writeFits (std::string const &fileName, std::string const &mode="w") const
 Write the object to a regular FITS file.
 
void writeFits (fits::MemFileManager &manager, std::string const &mode="w") const
 Write the object to a FITS image in memory.
 
void writeFits (fits::Fits &fitsfile) const
 Write the object to an already-open FITS object.
 

Static Public Member Functions

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

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.
 
void write (OutputArchiveHandle &handle) const override
 Write the object to one or more catalogs.
 
virtual std::string getPythonModule () const
 Return the fully-qualified Python module that should be imported to guarantee that its factory is registered.
 

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.
 

Related Symbols

(Note that these are not member symbols.)

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

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
table::Key< double > calibrationErr
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 }
bool isConstant() const
Is this photoCalib spatially constant?
Definition PhotoCalib.h:425

Member Function Documentation

◆ calibrateCatalog() [1/2]

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

Definition at line 350 of file PhotoCalib.cc.

350 {
351 std::vector<std::string> instFluxFields;
352 static std::string const SUFFIX = "_instFlux";
353 for (auto const &name : catalog.getSchema().getNames()) {
354 // Pick every field ending in "_instFlux", grabbing everything before that prefix.
355 if (name.size() > SUFFIX.size() + 1 &&
356 name.compare(name.size() - SUFFIX.size(), SUFFIX.size(), SUFFIX) == 0) {
357 instFluxFields.emplace_back(name.substr(0, name.size() - 9));
358 }
359 }
360 return calibrateCatalog(catalog, instFluxFields);
361}
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.
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 283 of file PhotoCalib.cc.

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

262 {
263 // Deep copy construct, as we're mutiplying in-place.
264 auto result = MaskedImage<float>(maskedImage, true);
265
266 if (_isConstant) {
267 *(result.getImage()) *= _calibrationMean;
268 } else {
269 _calibration->multiplyImage(*(result.getImage()), true); // only in the overlap region
270 }
271 if (includeScaleUncertainty) {
272 toNanojanskyVariance(maskedImage.getImage()->getArray(), maskedImage.getVariance()->getArray(),
273 _calibrationErr, result.getImage()->getArray(),
274 result.getVariance()->getArray());
275 } else {
276 toNanojanskyVariance(maskedImage.getImage()->getArray(), maskedImage.getVariance()->getArray(), 0,
277 result.getImage()->getArray(), result.getVariance()->getArray());
278 }
279
280 return result;
281}
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 239 of file PhotoCalib.cc.

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

◆ 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 222 of file PhotoCalib.cc.

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

◆ 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 226 of file PhotoCalib.cc.

226 {
227 throw LSST_EXCEPT(pex::exceptions::LogicError, "Not Implemented: See DM-10154.");
228}
#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.

◆ 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 253 of file PhotoCalib.cc.

253 {
254 return singleClassEquals(*this, other);
255}
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 439 of file PhotoCalib.h.

439{ 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 457 of file PhotoCalib.h.

457{ 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 510 of file PhotoCalib.cc.

510{ return getPhotoCalibPersistenceName(); }

◆ getPythonModule()

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

Return the fully-qualified Python module that should be imported to guarantee that its factory is registered.

Must be less than ArchiveIndexSchema::MAX_MODULE_LENGTH characters.

Will be ignored if empty.

Reimplemented in lsst::afw::image::FilterLabel, lsst::afw::cameraGeom::Detector, lsst::afw::cameraGeom::DetectorCollection, lsst::afw::cameraGeom::TransformMap, lsst::afw::detection::Footprint, lsst::afw::detection::GaussianPsf, lsst::afw::geom::SkyWcs, lsst::afw::geom::SpanSet, lsst::afw::geom::Transform< FromEndpoint, ToEndpoint >, lsst::afw::geom::Transform< afw::geom::Point2Endpoint, afw::geom::GenericEndpoint >, lsst::afw::image::ApCorrMap, lsst::afw::image::CoaddInputs, lsst::afw::image::TransmissionCurve, lsst::afw::math::ChebyshevBoundedField, lsst::afw::math::Function< double >, lsst::afw::math::Function< Kernel::Pixel >, lsst::afw::math::Function< Pixel >, lsst::afw::math::Function< ReturnT >, lsst::afw::math::Kernel, lsst::afw::math::PixelAreaBoundedField, lsst::afw::math::ProductBoundedField, lsst::afw::math::TransformBoundedField, lsst::afw::math::LanczosWarpingKernel, lsst::afw::math::BilinearWarpingKernel, lsst::afw::math::NearestWarpingKernel, lsst::afw::math::WarpingControl, lsst::afw::typehandling::StorableHelper< Base >, lsst::meas::algorithms::CoaddBoundedField, lsst::meas::algorithms::CoaddPsf, lsst::meas::algorithms::KernelPsf, lsst::meas::algorithms::WarpedPsf, lsst::meas::extensions::psfex::PsfexPsf, and lsst::meas::modelfit::Mixture.

Definition at line 36 of file Persistable.cc.

36{ return std::string(); }

◆ 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::geom::polygon::Polygon, lsst::afw::image::FilterLabel, lsst::afw::image::VisitInfo, and lsst::afw::typehandling::StorableHelper< Base >.

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 198 of file PhotoCalib.cc.

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

◆ 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 190 of file PhotoCalib.cc.

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

◆ 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 182 of file PhotoCalib.cc.

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

◆ instFluxToMagnitude() [4/7]

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

Definition at line 163 of file PhotoCalib.cc.

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

◆ instFluxToMagnitude() [5/7]

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

Definition at line 176 of file PhotoCalib.cc.

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

◆ 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 167 of file PhotoCalib.cc.

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

◆ 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 159 of file PhotoCalib.cc.

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

◆ 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 143 of file PhotoCalib.cc.

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

◆ 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 135 of file PhotoCalib.cc.

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

◆ 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 128 of file PhotoCalib.cc.

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

◆ instFluxToNanojansky() [4/7]

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

Definition at line 109 of file PhotoCalib.cc.

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

◆ instFluxToNanojansky() [5/7]

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

Definition at line 122 of file PhotoCalib.cc.

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

◆ 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 113 of file PhotoCalib.cc.

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

◆ 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 105 of file PhotoCalib.cc.

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

◆ isConstant()

bool lsst::afw::image::PhotoCalib::isConstant ( ) const
inline

Is this photoCalib spatially constant?

Returns
Is the calibration spatially constant?

Definition at line 425 of file PhotoCalib.h.

425{ return _isConstant; }

◆ 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 503 of file PhotoCalib.h.

503{ return true; }

◆ magnitudeToInstFlux() [1/2]

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

Definition at line 214 of file PhotoCalib.cc.

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

◆ 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 218 of file PhotoCalib.cc.

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

◆ 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 501 of file PhotoCalib.h.

501{ 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 230 of file PhotoCalib.cc.

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

◆ 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.

◆ 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.

◆ 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.

◆ 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 243 of file PhotoCalib.cc.

243 {
244 std::stringstream buffer;
245 if (_isConstant)
246 buffer << "spatially constant with ";
247 else
248 buffer << *_calibration << " with ";
249 buffer << "mean: " << _calibrationMean << " error: " << _calibrationErr;
250 return buffer.str();
251}
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 512 of file PhotoCalib.cc.

512 {
513 PhotoCalibSchema const &keys = PhotoCalibSchema::get();
514 table::BaseCatalog catalog = handle.makeCatalog(keys.schema);
515 auto record = catalog.addNew();
516 record->set(keys.calibrationMean, _calibrationMean);
517 record->set(keys.calibrationErr, _calibrationErr);
518 record->set(keys.isConstant, _isConstant);
519 record->set(keys.field, handle.put(_calibration));
520 record->set(keys.version, SERIALIZATION_VERSION);
521 handle.saveCatalog(catalog);
522}
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:489
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.

◆ 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 Symbol 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}

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