Namespace lsst::afw::geom

namespace geom

Typedefs

using TransformPoint2ToPoint2 = Transform<Point2Endpoint, Point2Endpoint>
using TransformPoint2ToGeneric = Transform<Point2Endpoint, GenericEndpoint>
using TransformPoint2ToSpherePoint = Transform<Point2Endpoint, SpherePointEndpoint>

Functions

std::ostream &operator<<(std::ostream &os, GenericEndpoint const &endpoint)

Print “GenericEndpoint(_n_)” to the ostream where _n_ is the number of axes, e.g. “GenericAxes(4)”.

std::ostream &operator<<(std::ostream &os, Point2Endpoint const &endpoint)

Print “Point2Endpoint()” to the ostream.

std::ostream &operator<<(std::ostream &os, SpherePointEndpoint const &endpoint)

Print “SpherePointEndpoint()” to the ostream.

Eigen::Matrix2d makeCdMatrix(lsst::geom::Angle const &scale, lsst::geom::Angle const &orientation = 0 * lsst::geom::degrees, bool flipX = false)

Make a WCS CD matrix

Return

the CD matrix, where element (i-1, j-1) corresponds to FITS keyword CDi_j and i, j have range [1, 2]

Parameters
  • [in] scale: Pixel scale as an angle on sky/pixels

  • [in] orientation: Position angle of pixel +Y, measured from N through E. At 0 degrees, +Y is along N and +X is along W/E if flipX false/true At 90 degrees, +Y is along E and +X is along N/S if flipX false/true

  • [in] flipX: Flip x axis? See orientation for details.

std::shared_ptr<SkyWcs> makeFlippedWcs(SkyWcs const &wcs, bool flipLR, bool flipTB, lsst::geom::Point2D const &center)

Return a copy of a FITS-WCS with pixel positions flipped around a specified center

Parameters
  • [in] wcs: The initial WCS

  • [in] flipLR: Flip pixel positions left/right about center

  • [in] flipTB: Flip pixel positions top/bottom about center

  • [in] center: Center pixel position of flip

std::shared_ptr<SkyWcs> makeModifiedWcs(TransformPoint2ToPoint2 const &pixelTransform, SkyWcs const &wcs, bool modifyActualPixels)

Create a new SkyWcs whose pixels are transformed by pixelTransform, as described below.

If modifyActualPixels is true and the cameraGeom::ACTUAL_PIXELS frame exists then pixelTransform is inserted just after the cameraGeom::ACTUAL_PIXELS frame:

newActualPixelsToPixels = pixelTransform -> oldActualPixelsToPixels

This is appropriate for shifting a WCS, e.g. when writing FITS metadata for a subimage.

If modifyActualPixels is false or the cameraGeom::ACTUAL_PIXELS frame does not exist then pixelTransform is inserted just after the cameraGeom::PIXELS frame:

newPixelsToIwc = pixelTransform -> oldPixelsToIwc

This is appropriate for inserting a model for optical distortion.

Other than the change described above, the new SkyWcs will be just like the old one.

Return

the new WCS

Parameters
  • [in] pixelTransform: Transform to insert

  • [in] wcs: Input WCS

  • [in] modifyActualPixels: Location at which to insert the transform; if true and the cameraGeom::ACTUAL_PIXELS frame is present then insert just after the cameraGeom::ACTUAL_PIXELS frame, else insert just after the cameraGeom::PIXELS frame.

std::shared_ptr<SkyWcs> makeSkyWcs(daf::base::PropertySet &metadata, bool strip = false)

Construct a SkyWcs from FITS keywords

This function is preferred over calling the SkyWcs metadata constructor directly because it allows us to change SkyWcs to an abstract base class in the future, without affecting code that constructs a WCS from FITS metadata.

Parameters
  • [in] metadata: FITS header metadata

  • [in] strip: If true: strip items from metadata used to create the WCS, such as RADESYS, EQUINOX, CTYPE12, CRPIX12, CRVAL12, etc. Always keep keywords that might be wanted for other purpposes, including NAXIS12 and date-related keywords such as “DATE-OBS” and “TIMESYS” (but not “EQUINOX”).

Exceptions
  • lsst::pex::exceptions::TypeError: if the metadata does not describe a celestial WCS.

std::shared_ptr<SkyWcs> makeSkyWcs(lsst::geom::Point2D const &crpix, lsst::geom::SpherePoint const &crval, Eigen::Matrix2d const &cdMatrix, std::string const &projection = "TAN")

Construct a simple FITS SkyWcs with no distortion

Parameters
  • [in] crpix: Center of projection on the CCD using the LSST convention: 0, 0 is the lower left pixel of the parent image

  • [in] crval: Center of projection on the sky

  • [in] cdMatrix: CD matrix, where element (i-1, j-1) corresponds to FITS keyword CDi_j and i, j have range [1, 2]. May be computed by calling makeCdMatrix.

  • [in] projection: The name of the FITS WCS projection, e.g. “TAN” or “STG”

std::shared_ptr<SkyWcs> makeSkyWcs(TransformPoint2ToPoint2 const &pixelsToFieldAngle, lsst::geom::Angle const &orientation, bool flipX, lsst::geom::SpherePoint const &boresight, std::string const &projection = "TAN")

Construct a FITS SkyWcs from camera geometry

Return

a SkyWcs whose sky origin is the boresight and pixel origin is focal plane (0, 0).

Note

Unlike makeCdMatrix, orientation is with respect to the focal plane axes, not the CCD axes. This is because field angle is defined with respect to focal plane axes.

Parameters
  • [in] pixelsToFieldAngle: Transformation from pixels to field angle (in radians).

  • [in] orientation: Position angle of focal plane +Y, measured from N through E at crval. At 0 degrees, +Y is along N and +X is along W/E if flipX false/true. At 90 degrees, +Y is along E and +X is along N/S if flipX false/true.

  • [in] flipX: Flip x axis? See orientation for details.

  • [in] boresight: ICRS sky position at the boresight (field angle (0, 0)).

  • [in] projection: The name of the FITS WCS projection, e.g. “TAN” or “STG”.

std::shared_ptr<SkyWcs> makeTanSipWcs(lsst::geom::Point2D const &crpix, lsst::geom::SpherePoint const &crval, Eigen::Matrix2d const &cdMatrix, Eigen::MatrixXd const &sipA, Eigen::MatrixXd const &sipB)

Construct a TAN-SIP SkyWcs with forward SIP distortion terms and an iterative inverse.

Parameters
  • [in] crpix: Center of projection on the CCD using the LSST convention: 0, 0 is the lower left pixel of the parent image

  • [in] crval: Center of projection on the sky

  • [in] cdMatrix: CD matrix, where element (i-1, j-1) corresponds to FITS keyword CDi_j and i, j have range [1, 2]. May be computed by calling makeCdMatrix.

  • [in] sipA: Forward distortion matrix for axis 1

  • [in] sipB: Forward distortion matrix for axis 2

std::shared_ptr<SkyWcs> makeTanSipWcs(lsst::geom::Point2D const &crpix, lsst::geom::SpherePoint const &crval, Eigen::Matrix2d const &cdMatrix, Eigen::MatrixXd const &sipA, Eigen::MatrixXd const &sipB, Eigen::MatrixXd const &sipAp, Eigen::MatrixXd const &sipBp)

Construct a TAN WCS with forward and inverse SIP distortion terms.

Parameters
  • [in] crpix: Center of projection on the CCD using the LSST convention: 0, 0 is the lower left pixel of the parent image

  • [in] crval: Center of projection on the sky

  • [in] cdMatrix: CD matrix, where element (i-1, j-1) corresponds to FITS keyword CDi_j and i, j have range [1, 2]. May be computed by calling makeCdMatrix.

  • [in] sipA: Forward distortion matrix for axis 1

  • [in] sipB: Forward distortion matrix for axis 2

  • [in] sipAp: Reverse distortion matrix for axis 1

  • [in] sipBp: Reverse distortion matrix for axis 2

std::shared_ptr<TransformPoint2ToPoint2> makeWcsPairTransform(SkyWcs const &src, SkyWcs const &dst)

A Transform obtained by putting two SkyWcs objects “back to back”.

Provides basic exception safety.

Return

a Transform whose forward transformation converts from src pixels to dst pixels, and whose inverse transformation converts in the opposite direction.

Parameters
  • src: the WCS for the source pixels

  • dst: the WCS for the destination pixels

std::shared_ptr<TransformPoint2ToSpherePoint> getIntermediateWorldCoordsToSky(SkyWcs const &wcs, bool simplify = true)

Return a transform from intermediate world coordinates to sky

std::shared_ptr<TransformPoint2ToPoint2> getPixelToIntermediateWorldCoords(SkyWcs const &wcs, bool simplify = true)

Return a transform from pixel coordinates to intermediate world coordinates

The pixel frame is is the base frame: cameraGeom::ACTUAL_PIXELS, if present, else cameraGeom::PIXELS.

std::ostream &operator<<(std::ostream &os, SkyWcs const &wcs)

Print a SkyWcs to an ostream (delegates to SkyWcs.toString()).

template<class FromEndpoint, class ToEndpoint>
std::ostream &operator<<(std::ostream &os, Transform<FromEndpoint, ToEndpoint> const &transform)

Print a Transform to an ostream

The format is “Transform<_fromEndpoint_, _toEndpoint_>” where fromEndpoint and toEndpoint are the appropriate endpoint printed to the ostream; for example “Transform<GenericEndpoint(4), Point2Endpoint()>”

lsst::geom::AffineTransform linearizeTransform(TransformPoint2ToPoint2 const &original, lsst::geom::Point2D const &inPoint)

Approximate a Transform by its local linearization.

Return

an lsst::geom::AffineTransform whose value and Jacobian at inPoint match those of original. It may be invertible; in general, linearizations are invertible if the Jacobian at inPoint is invertible.

Template Parameters
  • FromEndpointToEndpoint: The endpoints of the transform.

Parameters
  • original: the Transform to linearize

  • inPoint: the point at which a linear approximation is desired

Exceptions
  • pex::exceptions::InvalidParameterError: Thrown if original does not have a well-defined value and Jacobian at inPoint Not exception safe.

std::shared_ptr<TransformPoint2ToPoint2> makeTransform(lsst::geom::AffineTransform const &affine)

Wrap an lsst::geom::AffineTransform as a Transform.

Provides basic exception safety.

Return

a Transform that that maps any lsst::geom::Point2D x to affine(x). It shall be invertible iff affine is invertible.

Parameters
  • affine: The lsst::geom::AffineTransform to wrap.

std::shared_ptr<TransformPoint2ToPoint2> makeRadialTransform(std::vector<double> const &coeffs)

A purely radial polynomial distortion.

The Transform transforms an input \(x\) to

\[ \frac{x}{r} \sum_{i=1}^{N} \mathrm{coeffs[i]} \ r^i \]
where \(r\) is the magnitude of \(x\).

Return

the radial distortion represented by coeffs. The Transform shall have an inverse, which may be approximate.

Parameters
  • coeffs: radial polynomial coefficients. May be an empty vector to represent the identity transformation; otherwise must have size > 1, coeffs[0] = 0, and coeffs[1] 0.

Exceptions
  • pex::exceptions::InvalidParameterError: Thrown if coeffs does not have the required format. Provides basic exception safety.

std::shared_ptr<TransformPoint2ToPoint2> makeRadialTransform(std::vector<double> const &forwardCoeffs, std::vector<double> const &inverseCoeffs)

A purely radial polynomial distortion.

Similar to makeRadialTransform(std::vector<double> const &), but allows the user to provide an inverse.

Return

the radial distortion represented by coeffs. The Transform shall have an inverse, whose accuracy is determined by the relationship between forwardCoeffs and inverseCoeffs.

Parameters
  • forwardCoeffs: radial polynomial coefficients. May be an empty vector to represent the identity transformation; otherwise must have size > 1, coeffs[0] = 0, and coeffs[1] 0.

  • inverseCoeffs: coefficients for the inverse transform, as above. Does not need to have the same degree as forwardCoeffs, but either both must be empty or neither must be empty.

Exceptions
  • pex::exceptions::InvalidParameterError: Thrown if forwardCoeffs or inverseCoeffs does not have the required format. Provides basic exception safety.

std::shared_ptr<TransformPoint2ToPoint2> makeIdentityTransform()

Trivial Transform x x.

Provides basic exception safety.

Return

a Transform mapping any lsst::geom::Point2D to itself. The Transform’s inverse shall be itself.

std::shared_ptr<daf::base::PropertyList> createTrivialWcsMetadata(std::string const &wcsName, lsst::geom::Point2I const &xy0)
void deleteBasicWcsMetadata(daf::base::PropertySet &metadata, std::string const &wcsName)
Eigen::Matrix2d getCdMatrixFromMetadata(daf::base::PropertySet &metadata)

Read a CD matrix from FITS WCS metadata

The elements of the returned matrix are in degrees

Exceptions
  • pex::exceptions::TypeError: if no CD matrix coefficients found (missing coefficients are set to 0, as usual, but they cannot all be missing).

lsst::geom::Point2I getImageXY0FromMetadata(daf::base::PropertySet &metadata, std::string const &wcsName, bool strip = false)
Eigen::MatrixXd getSipMatrixFromMetadata(daf::base::PropertySet const &metadata, std::string const &name)
bool hasSipMatrix(daf::base::PropertySet const &metadata, std::string const &name)
std::shared_ptr<daf::base::PropertyList> makeSipMatrixMetadata(Eigen::MatrixXd const &matrix, std::string const &name)
std::shared_ptr<daf::base::PropertyList> makeSimpleWcsMetadata(lsst::geom::Point2D const &crpix, lsst::geom::SpherePoint const &crval, Eigen::Matrix2d const &cdMatrix, std::string const &projection = "TAN")

Make FITS metadata for a simple FITS WCS (one with no distortion).

This can also be used as a starting point for creating metadata for more sophisticated FITS WCS.

Parameters
  • [in] crpix: Center of projection on the CCD using the LSST convention: 0, 0 is the lower left pixel of the image

  • [in] crval: Center of projection on the sky

  • [in] cdMatrix: CD matrix where element (i-1, j-1) corresponds to FITS keyword CDi_j and i, j have range [1, 2]

  • [in] projection: The name of the FITS WCS projection, e.g. “TAN” or “STG”

std::shared_ptr<daf::base::PropertyList> makeTanSipMetadata(lsst::geom::Point2D const &crpix, lsst::geom::SpherePoint const &crval, Eigen::Matrix2d const &cdMatrix, Eigen::MatrixXd const &sipA, Eigen::MatrixXd const &sipB)

Make metadata for a TAN-SIP WCS without inverse matrices

Parameters
  • [in] crpix: Center of projection on the CCD using the LSST convention: 0, 0 is the lower left pixel of the parent image

  • [in] crval: Center of projection on the sky

  • [in] cdMatrix: CD matrix, where element (i-1, j-1) corresponds to FITS keyword CDi_j and i, j have range [1, 2]. May be computed by calling makeCdMatrix.

  • [in] sipA: Forward distortion matrix for axis 1

  • [in] sipB: Forward distortion matrix for axis 2

std::shared_ptr<daf::base::PropertyList> makeTanSipMetadata(lsst::geom::Point2D const &crpix, lsst::geom::SpherePoint const &crval, Eigen::Matrix2d const &cdMatrix, Eigen::MatrixXd const &sipA, Eigen::MatrixXd const &sipB, Eigen::MatrixXd const &sipAp, Eigen::MatrixXd const &sipBp)

Make metadata for a TAN-SIP WCS

Parameters
  • [in] crpix: Center of projection on the CCD using the LSST convention: 0, 0 is the lower left pixel of the parent image

  • [in] crval: Center of projection on the sky

  • [in] cdMatrix: CD matrix, where element (i-1, j-1) corresponds to FITS keyword CDi_j and i, j have range [1, 2]. May be computed by calling makeCdMatrix.

  • [in] sipA: Forward distortion matrix for axis 1

  • [in] sipB: Forward distortion matrix for axis 2

  • [in] sipAp: Reverse distortion matrix for axis 1

  • [in] sipBp: Reverse distortion matrix for axis 2

template<typename PointT, typename ArrayT>
class BaseEndpoint
#include <Endpoint.h>

Virtual base class for endpoints, which are helper classes for Transform

Endpoints transform points and lists of points from LSST-specific data types, such as lsst::geom::Point2D and lsst::geom::SpherePoint, to a form accepted by ast::Mapping.tran. Each type of endpoint is used for a particular LSST data type, for example:

Endpoints use the following forms of data for raw data:

  • std::vector<double> for a single point

  • ndarray<double, 2, 2> with dimensions number of axes x number of points for an array of points

Endpoints are designed as helper classes for Transform. Each transform has a two endpoints: one for input data and one for output data.

Endpoint also provides two methods to work with ast::Frames:

  • normalizeFrame verifies that a frame is the correct type, and adjusts its settings if necessary

  • makeFrame creates a new frame with the correct type and settings

Template Parameters
  • PointT: LSST data type for one point

  • ArrayT: LSST data type for an array of points

template<typename PointT>
class BaseVectorEndpoint : public lsst::afw::geom::BaseEndpoint<PointT, std::vector<PointT>>
#include <Endpoint.h>

Base class for endpoints with Array = std::vector<Point> where Point has 2 dimensions

Note

Subclasses must provide dataFromPoint, dataFromArray, pointFromData and arrayFromData

class GenericEndpoint : public lsst::afw::geom::BaseEndpoint<std::vector<double>, ndarray::Array<double, 2, 2>>
#include <Endpoint.h>

A generic endpoint for data in the format used by ast::Mapping

Thus supports all ast frame classes and any number of axes, and thus can be used as an endpoint for any ast::Mapping.

class Point2Endpoint : public lsst::afw::geom::BaseVectorEndpoint<lsst::geom::Point2D>
#include <Endpoint.h>

An endpoint for lsst::geom::Point2D

class SipApproximation
#include <SipApproximation.h>

A fitter and results class for approximating a general Transform in a form compatible with FITS WCS persistence.

The Simple Imaging Polynomial (SIP) convention (Shupe et al 2005) adds forward and reverse polynomial mappings to a standard projection FITS WCS projection (e.g. “TAN” for gnomonic) that relate Intermediate World Coordinates (see Calabretta & Greisen 2002) to image pixel coordinates. The SIP “forward” transform is defined by polynomial coeffients \(A\) and \(B\) that map pixel coordinates \((u, v)\) to Intermediate World Coordinates \((x, y)\) via

\[\begin{split} \boldsymbol{S}\left[\begin{array}{c} x \\ y \end{array}\right] \equiv \left[\begin{array}{c} x_s \\ y_s \end{array}\right] = \left[\begin{array}{c} (u - u_0) + \displaystyle\sum_{p,q}^{0 \le p + q \le N} \mathrm{A}_{p,q} (u - u_0)^p (v - v_0)^q \\ (v - v_0) + \displaystyle\sum_{p,q}^{0 \le p + q \le N} \mathrm{B}_{p,q} (u - u_0)^p (v - v_0)^q \end{array}\right] \end{split}\]
The reverse transform has essentially the same form:
\[\begin{split} \left[\begin{array}{c} u - u_0 \\ v - v_0 \end{array}\right] = \left[\begin{array}{c} x_s + \displaystyle\sum_{p,q}^{0 \le p + q \le N} \mathrm{AP}_{p,q} x_s^p y_s^q \\ y_s + \displaystyle\sum_{p,q}^{0 \le p + q \le N} \mathrm{BP}_{p,q} x_s^p y_s^q \end{array}\right] \end{split}\]
In both cases, \((u_0, v_0)\) is the pixel origin (CRPIX in FITS WCS) and \(\boldsymbol{S}\) is the inverse of the Jacobian “CD” matrix. Both CRPIX and CD are considered fixed inputs, and we do not attempt to null the zeroth- and first-order terms of \(A\) and \(B\) (as some SIP fitters do); together, these conventions make solving for the coefficients a much simpler linear problem.

While LSST WCSs are in general too complex to be described exactly in FITS WCS, they can generally be closely approximated by standard FITS WCS projection with additional SIP distortions. This class fits such an approximation, given a TransformPoint2ToPoint2 object that represents the exact mapping from pixels to Intermediate World Coordinates with a SIP distortion.

Note

In the implementation, we typically refer to \((u-u_0, v-v_0)\) as dpix (for “pixel delta”), and \((x_s, y_s)\) as siwc

(for “scaled

intermediate world coordinates”).

class SkyWcs : public lsst::afw::table::io::PersistableFacade<SkyWcs>, public Storable
#include <SkyWcs.h>

A 2-dimensional celestial WCS that transform pixels to ICRS RA/Dec, using the LSST standard for pixels.

SkyWcs is an immutable object that can not only represent any standard FITS WCS, but can also contain arbitrary Transforms, e.g. to model optical distortion or pixel imperfections.

In order to make a SkyWcs that models optical distortion, say, it is usually simplest to start with a standard FITS WCS (such as a TAN WCS) as an approximation, then insert a transform that models optical distortion by calling makeModifiedWcs. However, it is also possible to build a SkyWcs entirely from transforms, if you prefer, by building an ast::FrameDict and constructing the SkyWcs from that.

Frames of reference

SkyWcs internally keeps track of the following frames of reference:

  • cameraGeom::ACTUAL_PIXELS (optional): “actual” pixel position using the LSST standard. This has the same meaning as the cameraGeom::ACTUAL_PIXELS frame: actual pixels include effects such as “tree ring” distortions and electrical effects at the edge of CCDs. This frame should only be provided if there is a reasonable model for these imperfections.

  • cameraGeom::PIXELS: nominal pixel position, using the LSST standard. This has the same meaning as the cameraGeom::PIXELS frame: nominal pixels may be rectangular, but are uniform in size and spacing.

  • IWC: intermediate world coordinates (the FITS WCS concept).

  • SKY: position on the sky as ICRS, with standard RA, Dec axis order.

Pixel position standards

The LSST standard for pixel position is: 0,0 is the center of the lower left image pixel. The FITS standard for pixel position is: 1,1 is the center of the lower left image pixel.

LSST and FITS also use a different origin for subimages:

  • LSST pixel position is in the frame of reference of the parent image

  • FITS pixel position is in the frame of reference of the subimage However, SkyWcs does not keep track of this difference. Code that persists and unpersists SkyWcs using FITS-WCS header cards must handle the offset, e.g. by calling copyAtShiftedPixelOrigin

Internal details: the contained ast::FrameDict

SkyWcs contains an ast::FrameDict which transforms from pixels to sky (in radians) in the forward direction.

This FrameDict contains the named frames described in frames of reference, e.g. “SKY”, “IWC”, cameraGeom::PIXELS and possibly cameraGeom::ACTUAL_PIXELS. “SKY” is the current frame. If cameraGeom::ACTUAL_PIXELS is present then it is the base frame, otherwise cameraGeom::PIXELS is the base frame.

The “SKY” frame is of type ast::SkyFrame and has the following attributes:

  • SkyRef is set to the sky origin of the WCS (ICRS RA, Dec) in radians.

  • SkyRefIs is set to “Ignored” so that SkyRef is not used in transformations.

The other frames are of type ast::Frame and have 2 axes.

class Span
#include <Span.h>

A range of pixels within one row of an Image

class SpanPixelIterator : public boost::iterator_facade<SpanPixelIterator, lsst::geom::Point2I const, boost::random_access_traversal_tag>
#include <SpanPixelIterator.h>

An iterator that yields lsst::geom::Point2I and increases in the x direction.

This is used to iterate over the pixels in a Span, and by extension to iterate over regions like boxes and ellipses.

class SpherePointEndpoint : public lsst::afw::geom::BaseVectorEndpoint<lsst::geom::SpherePoint>
#include <Endpoint.h>

An endpoint for lsst::geom::SpherePoint

A SpherePointEndpoint always has 2 axes: longitude, latitude

template<class FromEndpoint, class ToEndpoint>
class Transform : public lsst::afw::table::io::PersistableFacade<Transform<FromEndpoint, ToEndpoint>>, public lsst::afw::table::io::Persistable
#include <Transform.h>

Transform LSST spatial data, such as lsst::geom::Point2D and lsst::geom::SpherePoint, using an AST mapping.

This class contains two Endpoints, to specify the “from” and “to” LSST data type, and an ast::Mapping to specify the transformation.

Depending on the ast::FrameSet or ast::Mapping used to define it, a Transform may provide either a forward transform, an inverse transform, or both. In particular, the inverse of a forward-only transform is an inverse-only transform. The hasForward and hasInverse methods can be used to check which transforms are available.

Unless otherwise stated, all constructors and methods may throw std::runtime_error to indicate internal errors within AST.

Transforms are always immutable.

Note

You gain some safety by constructing a Transform from an ast::FrameSet, since the base and current frames in the FrameSet can be checked against by the appropriate endpoint.

Note

”In place” versions of applyForward and applyInverse are not available because data must be copied when converting from LSST data types to the type used by astshim, so it didn’t seem worth the bother.

namespace detail

Functions

std::shared_ptr<ast::FrameSet> readFitsWcs(daf::base::PropertySet &metadata, bool strip = true)

Read a FITS convention WCS FrameSet from FITS metadata

The resulting FrameSet may be any kind of WCS supported by FITS; if it is a celestial WCS then 1,1 will be the lower left corner of the image (the FITS convention, not the LSST convention).

This routine replaces RADECSYS with RADESYS if the former is present and the latter is not, since that is a common misspelling in FITS headers.

The returned FrameSet will have an IWC (intermediate world coordinate system) frame.

Parameters
  • [inout] metadata: FITS header cards

  • [in] strip: If true: strip items from metadata used to create the WCS, such as RADESYS, EQUINOX, CTYPE12, CRPIX12, CRVAL12, etc. Always keep keywords that might be wanted for other purpposes, including NAXIS12 and date-related keywords such as “DATE-OBS” and “TIMESYS” (but not “EQUINOX”).

Exceptions
  • pex::exceptions::TypeError: if the metadata does not contain a FITS-WCS

std::shared_ptr<ast::FrameDict> readLsstSkyWcs(daf::base::PropertySet &metadata, bool strip = true)

Read an LSST celestial WCS FrameDict from a FITS header.

Calls getImageXY0FromMetadata to determine image XY0.

Saves CRVAL by setting the output SkyFrame’s SkyRef to CRVAL and SkyRefIs=”Ignore” (so SkyRef is not treated as an offset).

The frames of the returned WCS will be as follows:

  • ”PIXELS” (base frame): pixel frame with 0,0 the lower left corner of the image (LSST convention)

  • ”IWC”: FITS WCS intermediate world coordinate system

  • ”SKY” (current frame): an ast::SkyFrame with domain “SKY”: ICRS RA, Dec

Warning

Does not compensate for the offset between a subimage and its parent image; callers must do that manually, e.g. by calling SkyWcs::copyAtShiftedPosition.

Warning

the general SkyWcs generated by LSST software cannot be exactly represented using standard WCS FITS cards, and so this function cannot read those. This function is intended for two purposes:

  • Read the standard FITS WCS found in raw data and other images from non-LSST observatories and convert it to the LSST pixel convention.

  • Read the approximate FITS WCS that LSST writes to FITS images (for use by non-LSST code).

All frames are instances of ast::Frame except the SKY frame. All have 2 axes.

Parameters
  • [inout] metadata: FITS header cards

  • [in] strip: If true: strip items from metadata used to create the WCS, such as RADESYS, EQUINOX, CTYPE12, CRPIX12, CRVAL12, etc. Always keep keywords that might be wanted for other purpposes, including NAXIS12 and date-related keywords such as “DATE-OBS” and “TIMESYS” (but not “EQUINOX”).

Exceptions
  • lsst::pex::exceptions::TypeError: if the metadata does not describe a celestial WCS.

std::shared_ptr<daf::base::PropertyList> getPropertyListFromFitsChan(ast::FitsChan &fitsChan)

Copy values from an AST FitsChan into a PropertyList

Warning

COMMENT and HISTORY cards are treated as string values

Exceptions
  • lsst::pex::exceptions::TypeError: if fitsChan contains cards whose type is not supported by PropertyList: complex numbers, or cards with no value

ast::FitsChan getFitsChanFromPropertyList(daf::base::PropertySet &metadata, std::set<std::string> const &excludeNames = {}, std::string options = "")

Construct AST FitsChan from PropertyList

Return

an ast::FitsChan representing this PropertyList

Parameters
  • [in] metadata: Metadata to transfer; if this is a PropertyList then the order of items is preserved.

  • [in] excludeNames: Names of entries to exclude from the returned header string.

  • [in] options: Options string to pass to ast::FitsChan constructor.

template<class Transform>
std::shared_ptr<Transform> readStream(std::istream &is)

Deserialize a Transform from an input stream

Template Parameters
  • Transform: the Transform class; can be Transform<FromEndpoint, ToEndpoint>, SkyWcs, or any other compatible class, i.e. it must support the following (see Transform.h for details):

    • a constructor that takes an ast::FrameSet

    • static method getShortClassName

    • method getMapping

Parameters
  • [in] is: input stream from which to deserialize this Transform

template<class Transform>
void writeStream(Transform const &transform, std::ostream &os)

Serialize a Transform to an output stream

Version 1 format is as follows:

  • The version number (an integer)

  • A space

  • The short class name, as obtained from getShortClassName

  • A space

  • The contained ast::FrameSet written using FrameSet.show(os, false)

Parameters
  • [out] os: output stream to which to serialize this Transform

  • [in] transform: Transform to serialize

Variables

constexpr int serializationVersion = 1

version of serialization used when writing (older versions may also be supported when reading)

namespace ellipses

Typedefs

typedef Separable<Distortion, DeterminantRadius> SeparableDistortionDeterminantRadius
typedef Separable<Distortion, TraceRadius> SeparableDistortionTraceRadius
typedef Separable<Distortion, LogDeterminantRadius> SeparableDistortionLogDeterminantRadius
typedef Separable<Distortion, LogTraceRadius> SeparableDistortionLogTraceRadius
typedef Separable<ConformalShear, DeterminantRadius> SeparableConformalShearDeterminantRadius
typedef Separable<ConformalShear, TraceRadius> SeparableConformalShearTraceRadius
typedef Separable<ConformalShear, LogDeterminantRadius> SeparableConformalShearLogDeterminantRadius
typedef Separable<ConformalShear, LogTraceRadius> SeparableConformalShearLogTraceRadius
typedef Separable<ReducedShear, DeterminantRadius> SeparableReducedShearDeterminantRadius
typedef Separable<ReducedShear, TraceRadius> SeparableReducedShearTraceRadius
typedef Separable<ReducedShear, LogDeterminantRadius> SeparableReducedShearLogDeterminantRadius
typedef Separable<ReducedShear, LogTraceRadius> SeparableReducedShearLogTraceRadius
class Axes : public lsst::afw::geom::ellipses::BaseCore
#include <Axes.h>

An ellipse core for the semimajor/semiminor axis and position angle parametrization (a,b,theta).

class BaseCore
#include <BaseCore.h>

A base class for parametrizations of the “core” of an ellipse - the ellipticity and size.

A subclass of BaseCore provides a particular interpretation of the three pointing point values that define an ellipse’s size and ellipticity (including position angle). All core subclasses are implicitly convertible and can be assigned to from any other core.

Subclassed by lsst::afw::geom::ellipses::Axes, lsst::afw::geom::ellipses::Quadrupole, lsst::afw::geom::ellipses::Separable< Ellipticity_, Radius_ >

Coordinate transforms

These member functions transform the ellipse by the given lsst::geom::LinearTransform. The transform can be done in-place by calling inPlace() on the returned expression object, or returned as a new shared_ptr by calling copy().

class ConformalShear : public lsst::afw::geom::ellipses::detail::EllipticityBase
#include <ConformalShear.h>

A logarithmic complex ellipticity with magnitude \(|e| = \ln (a/b) \).

For a more complete definition, see Bernstein and Jarvis (2002); this the same as their conformal shear \(\eta\) (eq. 2.3-2.6).

class DeterminantRadius
#include <radii.h>

The radius defined as the 4th root of the determinant of the quadrupole matrix.

The determinant radius is equal to the standard radius for a circle, and \(\pi R_{det}^2\) is the area of the ellipse.

class Distortion : public lsst::afw::geom::ellipses::detail::EllipticityBase
#include <Distortion.h>

A complex ellipticity with magnitude \(|e| = \frac{a^2 - b^2}{a^2 + b^2}\).

For a more complete definition, see Bernstein and Jarvis (2002); this the same as their distortion \(\delta\) (eq. 2.7).

class Ellipse
#include <Ellipse.h>

An ellipse defined by an arbitrary BaseCore and a center point.

An ellipse is composed of its center coordinate and its Core - a parametrization of the ellipticity and size of the ellipse. Setting the core of an ellipse never changes the type of the contained core, it merely sets the parameters of that core by converting the parameters.

Coordinate transforms

These member functions transform the ellipse by the given lsst::geom::AffineTransform. The transform can be done in-place by calling inPlace() on the returned expression object, or returned as a new shared_ptr by calling copy().

class LogDeterminantRadius
#include <radii.h>

The natural logarithm of the DeterminantRadius

class LogTraceRadius
#include <radii.h>

The natural logarithm of the TraceRadius

class Parametric
#include <Parametric.h>

A functor that returns points on the boundary of the ellipse as a function of a parameter that runs between 0 and 2 pi (but is not angle).

class PixelRegion
#include <PixelRegion.h>

A pixelized region containing all pixels whose centers are within an Ellipse.

The pixel region for an ellipse may be larger or smaller in area than the ellipse itself, depending on the details of where pixel centers land, and it may be empty even if the area of the ellipse is nonzero.

class Quadrupole : public lsst::afw::geom::ellipses::BaseCore
#include <Quadrupole.h>

An ellipse core with quadrupole moments as parameters.

class ReducedShear : public lsst::afw::geom::ellipses::detail::EllipticityBase
#include <ReducedShear.h>

A complex ellipticity with magnitude \(|e| = \frac{a-b}{a+b} \).

For a more complete definition, see Bernstein and Jarvis (2002); this the same as their reduced shear \(g\) (eq. 2.8).

template<typename Ellipticity_, typename Radius_>
class Separable : public lsst::afw::geom::ellipses::BaseCore
#include <Separable.h>

An ellipse core with a complex ellipticity and radius parameterization.

class TraceRadius
#include <radii.h>

The radius defined as \(\sqrt{0.5(I_{xx} + I_{yy})}\)

The trace radius is equal to the standard radius for a circle

namespace detail
class EllipticityBase
#include <EllipticityBase.h>

EllipticityBase is a base class for complex ellipticity types.

EllipticityBase does not have a virtual destructor, and only exists for code reuse purposes. The ellipticity classes are not polymorphic simply to keep them small.

Subclassed by lsst::afw::geom::ellipses::ConformalShear, lsst::afw::geom::ellipses::Distortion, lsst::afw::geom::ellipses::ReducedShear

namespace polygon

Functions

LSST_EXCEPTION_TYPE(SinglePolygonException, lsst::pex::exceptions::RuntimeError, lsst::afw::geom::polygon::SinglePolygonException)

An exception that indicates the single-polygon assumption has been violated

The single-polygon assumption is used in Polygon::intersectionSingle and Polygon::unionSingle.

std::ostream &operator<<(std::ostream &os, Polygon const &poly)

Stream polygon.

class Polygon : public lsst::afw::table::io::PersistableFacade<Polygon>, public Storable
#include <Polygon.h>

Cartesian polygons

Polygons are defined by a set of vertices