File PhotometryTransform.h

namespace lsst

Class for a simple mapping implementing a generic AstrometryTransform.

Remove all non-astronomical counts from the Chunk Exposure’s pixels.

Forward declarations for lsst::utils::Cache

For details on the Cache class, see the Cache.h file.

It uses a template rather than a pointer so that the derived classes can use the specifics of the transform. The class simplePolyMapping overloads a few routines.

A base class for image defects

Numeric constants used by the Integrate.h integrator routines.

Compute Image Statistics

Note

Gauss-Kronrod-Patterson quadrature coefficients for use in quadpack routine qng. These coefficients were calculated with 101 decimal digit arithmetic by L. W. Fullerton, Bell Labs, Nov 1981.

Note

The Statistics class itself can only handle lsst::afw::image::MaskedImage() types. The philosophy has been to handle other types by making them look like lsst::afw::image::MaskedImage() and reusing that code. Users should have no need to instantiate a Statistics object directly, but should use the overloaded makeStatistics() factory functions.

namespace jointcal
class FluxTransformChebyshev : public lsst::jointcal::PhotometryTransformChebyshev
#include <PhotometryTransform.h>

nth-order 2d Chebyshev photometry transform, times the input flux.

Public Functions

FluxTransformChebyshev(size_t order, geom::Box2D const &bbox)
FluxTransformChebyshev(ndarray::Array<double, 2, 2> const &coefficients, geom::Box2D const &bbox)
double transform(double x, double y, double value) const

Return the transform of value at (x,y).

void computeParameterDerivatives(double x, double y, double value, Eigen::Ref<Eigen::VectorXd> derivatives) const

Compute the derivatives with respect to the parameters (i.e. the coefficients).

Parameters
  • [in] x: The x coordinate to compute at (in the appropriate units for this transform).

  • [in] y: The y coordinate to compute at (in the appropriate units for this transform).

  • [in] value: The instrument flux or magnitude to compute the derivative at.

  • [out] derivatives: The computed derivatives, in the same order as the deltas in offsetParams.

std::shared_ptr<PhotometryTransform> clone() const

return a copy (allocated by new) of the transformation.

class FluxTransformSpatiallyInvariant : public lsst::jointcal::PhotometryTransformSpatiallyInvariant
#include <PhotometryTransform.h>

Photometric offset independent of position, defined as (fluxMag0)^-1.

initialCalibFlux * SpatiallyInvariantTransform -> correctedFlux

Public Functions

FluxTransformSpatiallyInvariant(double value = 1)
double transform(double x, double y, double value) const

Return the transform of value at (x,y).

double transformError(double x, double y, double value, double valueErr) const

Return the transformed valueErr at Point(x,y).

std::shared_ptr<PhotometryTransform> clone() const

return a copy (allocated by new) of the transformation.

void computeParameterDerivatives(double x, double y, double value, Eigen::Ref<Eigen::VectorXd> derivatives) const

Compute the derivatives with respect to the parameters (i.e. the coefficients).

Parameters
  • [in] x: The x coordinate to compute at (in the appropriate units for this transform).

  • [in] y: The y coordinate to compute at (in the appropriate units for this transform).

  • [in] value: The instrument flux or magnitude to compute the derivative at.

  • [out] derivatives: The computed derivatives, in the same order as the deltas in offsetParams.

class MagnitudeTransformChebyshev : public lsst::jointcal::PhotometryTransformChebyshev
#include <PhotometryTransform.h>

nth-order 2d Chebyshev photometry transform, plus the input flux.

Public Functions

MagnitudeTransformChebyshev(size_t order, geom::Box2D const &bbox)
MagnitudeTransformChebyshev(ndarray::Array<double, 2, 2> const &coefficients, geom::Box2D const &bbox)
double transform(double x, double y, double value) const

Return the transform of value at (x,y).

void computeParameterDerivatives(double x, double y, double value, Eigen::Ref<Eigen::VectorXd> derivatives) const

Compute the derivatives with respect to the parameters (i.e. the coefficients).

Parameters
  • [in] x: The x coordinate to compute at (in the appropriate units for this transform).

  • [in] y: The y coordinate to compute at (in the appropriate units for this transform).

  • [in] value: The instrument flux or magnitude to compute the derivative at.

  • [out] derivatives: The computed derivatives, in the same order as the deltas in offsetParams.

std::shared_ptr<PhotometryTransform> clone() const

return a copy (allocated by new) of the transformation.

class MagnitudeTransformSpatiallyInvariant : public lsst::jointcal::PhotometryTransformSpatiallyInvariant
#include <PhotometryTransform.h>

Photometric offset independent of position, defined as -2.5 * log(flux / fluxMag0).

initialMagnitude + SpatiallyInvariantTransform -> correctedMagnitude

Public Functions

MagnitudeTransformSpatiallyInvariant(double value = 0)
double transform(double x, double y, double mag) const

Return the transform of value at (x,y).

double transformError(double x, double y, double value, double valueErr) const

Return the transformed valueErr at Point(x,y).

std::shared_ptr<PhotometryTransform> clone() const

return a copy (allocated by new) of the transformation.

void computeParameterDerivatives(double x, double y, double value, Eigen::Ref<Eigen::VectorXd> derivatives) const

Compute the derivatives with respect to the parameters (i.e. the coefficients).

Parameters
  • [in] x: The x coordinate to compute at (in the appropriate units for this transform).

  • [in] y: The y coordinate to compute at (in the appropriate units for this transform).

  • [in] value: The instrument flux or magnitude to compute the derivative at.

  • [out] derivatives: The computed derivatives, in the same order as the deltas in offsetParams.

class PhotometryTransform
#include <PhotometryTransform.h>

A photometric transform, defined in terms of the input flux or magnitude.

Unit agnostic: a higher level Model must keep track of the units going into and out of of its Transforms.

See

lsst::afw::image::PhotoCalib

Subclassed by lsst::jointcal::PhotometryTransformChebyshev, lsst::jointcal::PhotometryTransformSpatiallyInvariant

Public Functions

virtual double transform(double x, double y, double value) const = 0

Return the transform of value at (x,y).

double transform(Point const &in, double value) const

Return the transformed value at Point(x,y).

virtual double transformError(double x, double y, double value, double valueErr) const = 0

Return the transformed valueErr at Point(x,y).

double transformError(Point const &in, double value, double valueErr) const

Return the transformed valueErr at Point(x,y).

virtual void dump(std::ostream &stream = std::cout) const = 0

dumps the transform coefficients to stream.

virtual std::size_t getNpar() const = 0

Return the number of parameters (used to compute chisq)

virtual void offsetParams(Eigen::VectorXd const &delta) = 0

Offset the parameters by some (negative) amount during fitting.

Equivalent to flatten(parameters) -= delta

Ordering of delta is the same as the ordering of the derivatives returned from computeParameterDerivatives.

virtual std::shared_ptr<PhotometryTransform> clone() const = 0

return a copy (allocated by new) of the transformation.

virtual void computeParameterDerivatives(double x, double y, double value, Eigen::Ref<Eigen::VectorXd> derivatives) const = 0

Compute the derivatives with respect to the parameters (i.e. the coefficients).

Parameters
  • [in] x: The x coordinate to compute at (in the appropriate units for this transform).

  • [in] y: The y coordinate to compute at (in the appropriate units for this transform).

  • [in] value: The instrument flux or magnitude to compute the derivative at.

  • [out] derivatives: The computed derivatives, in the same order as the deltas in offsetParams.

virtual Eigen::VectorXd getParameters() const = 0

Get a copy of the parameters of this model, in the same order as offsetParams.

Friends

std::ostream &operator<<(std::ostream &s, PhotometryTransform const &transform)
class PhotometryTransformChebyshev : public lsst::jointcal::PhotometryTransform
#include <PhotometryTransform.h>

nth-order 2d Chebyshev photometry transform.

The 2-d Chebyshev polynomial used here is defined as:

\[ f(x,y) = \sum_i \sum_j a_{i,j} T_i(x) T_j(y) \]

where \(T_n(x)\) is the n-th order Chebyshev polynomial of \(x\) and \(a_{i,j}\) is the corresponding coefficient of the (i,j) polynomial term.

Note that the polynomial order=n means that the highest terms will be of the form:

\[ a_{0,n}*x^n*y^0, a_{n-1,1}*x^(n-1)*y^1, ..., a_{1,n-1}*x^1*y^(n-1), a_{n,0}*x^0*y^n \]

Subclassed by lsst::jointcal::FluxTransformChebyshev, lsst::jointcal::MagnitudeTransformChebyshev

Public Functions

PhotometryTransformChebyshev(size_t order, geom::Box2D const &bbox, bool identity)

Create a Chebyshev transform with terms up to order in (x*y).

Parameters
  • [in] order: The maximum order in (x*y).

  • [in] bbox: The bounding box it is valid within, to rescale it to [-1,1].

  • [in] identity: If true, set a_0,0==1, otherwise all coefficients are 0.

PhotometryTransformChebyshev(ndarray::Array<double, 2, 2> const &coefficients, geom::Box2D const &bbox)

Create a Chebyshev transform with the specified coefficients.

The polynomial order is determined from the number of coefficients, taking only the anti-diagonal upper triangle portion of the passed-in coefficients

Parameters
  • coefficients: The polynomial coefficients.

  • [in] bbox: The bounding box it is valid within, to rescale it to [-1,1].

double transformError(double x, double y, double value, double valueErr) const

Return the transformed valueErr at Point(x,y).

void dump(std::ostream &stream = std::cout) const

dumps the transform coefficients to stream.

std::size_t getNpar() const

Return the number of parameters (used to compute chisq)

void offsetParams(Eigen::VectorXd const &delta)

Offset the parameters by some (negative) amount during fitting.

Equivalent to flatten(parameters) -= delta

Ordering of delta is the same as the ordering of the derivatives returned from computeParameterDerivatives.

ndarray::Array<double, 2, 2> getCoefficients() const

Get a copy of the coefficients of the polynomials, as a 2d array (NOTE: layout is [y][x])

Eigen::VectorXd getParameters() const

Get a copy of the parameters of this model, in the same order as offsetParams.

ndarray::Size getOrder() const
geom::Box2D getBBox() const
double mean(geom::Box2D const &bbox) const

Compute the mean of this tranform on the bbox (default to our bbox).

double mean() const

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

double integrate(geom::Box2D const &bbox) const
double integrate() const

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

Protected Functions

double computeChebyshev(double x, double y) const

Return the value of this polynomial at x,y. For use in the sublcass transform() methods.

void computeChebyshevDerivatives(double x, double y, Eigen::Ref<Eigen::VectorXd> derivatives) const

Set the derivatives of this polynomial at x,y. For use in the sublcass computeParameterDerivatives() methods.

Private Functions

double oneIntegral(double x, double y) const

Evaluate one limit of a definite 2-d integral (sum four of these to get the full integral).

Private Members

geom::Box2D _bbox
geom::AffineTransform _toChebyshevRange
ndarray::Array<double, 2, 2> _coefficients
ndarray::Size _order
ndarray::Size _nParameters
class PhotometryTransformSpatiallyInvariant : public lsst::jointcal::PhotometryTransform
#include <PhotometryTransform.h>

Photometry offset independent of position. Abstract class.

Subclassed by lsst::jointcal::FluxTransformSpatiallyInvariant, lsst::jointcal::MagnitudeTransformSpatiallyInvariant

Public Functions

PhotometryTransformSpatiallyInvariant(double value)
void dump(std::ostream &stream = std::cout) const

dumps the transform coefficients to stream.

std::size_t getNpar() const

Return the number of parameters (used to compute chisq)

void offsetParams(Eigen::VectorXd const &delta)

Offset the parameters by some (negative) amount during fitting.

Equivalent to flatten(parameters) -= delta

Ordering of delta is the same as the ordering of the derivatives returned from computeParameterDerivatives.

Eigen::VectorXd getParameters() const

Get a copy of the parameters of this model, in the same order as offsetParams.

Protected Functions

double getValue() const

Private Members

double _value

value of this transform at all locations.