Namespace lsst::jointcal

namespace jointcal

Typedefs

using RefFluxMapType = std::map<std::string, std::vector<double>>
typedef void() lsst::jointcal::AstrometryTransformFun(const double, const double, double &, double &, const void *)

signature of the user-provided routine that actually does the coordinate transform for UserTransform.

typedef StarList<BaseStar> BaseStarList
typedef BaseStarList::const_iterator BaseStarCIterator
typedef BaseStarList::iterator BaseStarIterator
typedef std::list<std::shared_ptr<CcdImage>> CcdImageList
typedef int VisitIdType
typedef int CcdIdType
typedef FittedStarList::const_iterator FittedStarCIterator
typedef FittedStarList::iterator FittedStarIterator
typedef MeasuredStarList::const_iterator MeasuredStarCIterator
typedef MeasuredStarList::iterator MeasuredStarIterator
typedef RefStarList::const_iterator RefStarCIterator
typedef RefStarList::iterator RefStarIterator
typedef std::list<StarMatch>::iterator StarMatchIterator
typedef std::list<StarMatch>::const_iterator StarMatchCIterator
typedef Eigen::Triplet<double> Trip

Enums

enum MinimizeResult

Return value of minimize()

Values:

Converged
Chi2Increased
Failed
NonFinite

Functions

std::ostream &operator<<(std::ostream &stream, AstrometryTransform const &transform)

Delegates to transform.dump()

std::unique_ptr<AstrometryTransform> compose(AstrometryTransform const &left, AstrometryTransform const &right)

Returns a pointer to a composition of transforms, representing left(right()).

Deletion of returned value to be done by caller.

If left->composeAndReduce(right) returns NULL, build a AstrometryTransformComposition and return it. This routine implements “run-time” compositions. When there is a possible “reduction” (e.g. compositions of polynomials), compose detects it and returns a genuine AstrometryTransform.

Return

The composed transform.

std::unique_ptr<AstrometryTransform> compose(AstrometryTransform const &left, AstrometryTransformIdentity const &right)
bool isIntegerShift(const AstrometryTransform *transform)

Shorthand test to tell if a transform is a simple integer shift.

std::shared_ptr<AstrometryTransformPolynomial> inversePolyTransform(AstrometryTransform const &forward, Frame const &domain, double const precision, std::size_t maxOrder = 9, std::size_t nSteps = 50)

Approximate the inverse by a polynomial, to some precision.

Return

A polynomial that best approximates forward.

Parameters
  • forward: Transform to be inverted.

  • [in] domain: The domain of forward.

  • [in] precision: Require that \(chi2 / (nsteps^2) < precision^2\).

  • [in] maxOrder: The maximum order allowed of the inverse polynomial.

  • [in] nSteps: The number of sample points per axis (nSteps^2 total points).

AstrometryTransformLinear normalizeCoordinatesTransform(const Frame &frame)
std::unique_ptr<AstrometryTransform> astrometryTransformRead(const std::string &fileName)

The virtual constructor from a file.

std::unique_ptr<AstrometryTransform> astrometryTransformRead(std::istream &s)

The virtual constructor from a file.

std::ostream &operator<<(std::ostream &out, CcdImageKey const &key)
BaseStarList &Fitted2Base(FittedStarList &This)
BaseStarList *Fitted2Base(FittedStarList *This)
const BaseStarList &Fitted2Base(const FittedStarList &This)
const BaseStarList *Fitted2Base(const FittedStarList *This)
std::unique_ptr<StarMatchList> matchSearchRotShift(BaseStarList &list1, BaseStarList &list2, const MatchConditions &conditions)

searches a geometrical transformation that goes from list1 to list2.

The found transformation is a field of the returned object, as well as the star pairs (the matches) that were constructed. (see StarMatchList class definition for more details). The various cuts are contained in conditions (see listmatch.h) for its contents. This routine searches a transformation that involves a shift and a rotation.

std::unique_ptr<StarMatchList> matchSearchRotShiftFlip(BaseStarList &list1, BaseStarList &list2, const MatchConditions &conditions)

same as above but searches also a flipped solution.

std::unique_ptr<StarMatchList> listMatchCollect(const BaseStarList &list1, const BaseStarList &list2, const AstrometryTransform *guess, const double maxDist)

assembles star matches.

It picks stars in list1, transforms them through guess, and collects closest star in list2, and builds a match if closer than maxDist).

std::unique_ptr<StarMatchList> listMatchCollect(const BaseStarList &list1, const BaseStarList &list2, const double maxDist)

same as before except that the transform is the identity

std::unique_ptr<AstrometryTransformLinear> listMatchupShift(const BaseStarList &list1, const BaseStarList &list2, const AstrometryTransform &transform, double maxShift, double binSize = 0)

searches for a 2 dimensional shift using a very crude histogram method.

std::unique_ptr<AstrometryTransform> listMatchCombinatorial(const BaseStarList &list1, const BaseStarList &list2, const MatchConditions &conditions = MatchConditions())
std::unique_ptr<AstrometryTransform> listMatchRefine(const BaseStarList &list1, const BaseStarList &list2, std::unique_ptr<AstrometryTransform> transform, const int maxOrder = 3)
BaseStarList &Measured2Base(MeasuredStarList &This)
BaseStarList *Measured2Base(MeasuredStarList *This)
const BaseStarList &Measured2Base(const MeasuredStarList &This)
const BaseStarList *Measured2Base(const MeasuredStarList *This)
std::ostream &operator<<(std::ostream &stream, const Point &point)

BaseStarList &Ref2Base(RefStarList &This)
BaseStarList *Ref2Base(RefStarList *This)
const BaseStarList &Ref2Base(const RefStarList &This)
const BaseStarList *Ref2Base(const RefStarList *This)
template<class Star>
std::ostream &operator<<(std::ostream &stream, const StarList<Star> &list)

enables

std::cout << my_list;

bool compareStar1(const StarMatch &one, const StarMatch &two)
bool sameStar1(const StarMatch &one, const StarMatch &two)
bool compareStar2(const StarMatch &one, const StarMatch &two)
bool sameStar2(const StarMatch &one, const StarMatch &two)
std::ostream &operator<<(std::ostream &stream, const StarMatch &match)
std::ostream &operator<<(std::ostream &stream, const StarMatchList &starMatchList)

A std::list of star matches,.

To be used as the argument to AstrometryTransform::fit routines. There is as well a StarMatch::fit routine which fits a polynomial by default, although the transform may be user-provided. The StarMatchList::refineTransform is a convenient tool to reject outliers. Given two catalogs, one can assemble a StarMatchList using utilities such as listMatchCollect. StarMatchList’s have write capabilities. NStarMatchList is a generalization of this 2-match to n-matches.

double computeDist2(const StarMatchList &S, const AstrometryTransform &transform)

sum of distance squared

double computeChi2(const StarMatchList &L, const AstrometryTransform &transform)

the actual chi2

class Associations
#include <Associations.h>

The class that implements the relations between MeasuredStar and FittedStar.

class AstrometryFit : public lsst::jointcal::FitterBase
#include <AstrometryFit.h>

Class that handles the astrometric least squares problem.

This is the class that actually computes the quantities required to carry out a LS astrometric fit wrt distortion mappings and coordinates of common objects. Namely it computes the Jacobian and gradient of the chi2 (w.r.t. parameters), and the Chi2 itself. It interfaces with the actual modelling of distortions via a mimimum virtual interface AstrometryModel, and the actual mappings via an other virtual interface : Mapping.

In short AstrometryFit aims at computing derivatives of least quares. The terms of the chi2 are of two kinds:

kind 1 -> (T(X_M) - p(F))^T W (T(X_M) - p(F))

with X_M is a measured (2d) position in CCD coordinates, F refers to the position of the object in some space, defined in practise by p. There is one such term per measurement. The default setup would be that p is the projection from sky to some tangent plane and hence T maps the CCD coordinates onto this TP. p is obtained via the DistorsionModel and can be different for all CcdImage’s. Depending on what is beeing fitted, one could imagine cases where the projector p is the same for all CcdImages.

Kind 2 -> (p’(F)-p’(R))^T W_R (p’(F)-p’(R)) R refers to some externally-provided reference object position, and p’ to some projector from sky to some plane. The reference objects define the overall coordinate frame, which is required when all T and all F are fitted simultaneously. There is one such term per external reference object. There can be more F (fitted) objects than R (reference) objects.

In the same framework, one can fit relative transforms between images by setting p = Identity for all input CcdImages and not fitting T for one of the CcdImage’s. One does not need reference object and would then naturally not have any Kind 2 terms.

class AstrometryMapping
#include <AstrometryMapping.h>

virtual class needed in the abstraction of the distortion model

Subclassed by lsst::jointcal::ChipVisitAstrometryMapping, lsst::jointcal::SimpleAstrometryMapping

class AstrometryModel
#include <AstrometryModel.h>

Interface between AstrometryFit and the combinations of Mappings from pixels to some tangent plane (aka distortions).

Parameters
  • log: Logger to send messages to, to keep names consistent when logging.

Subclassed by lsst::jointcal::ConstrainedAstrometryModel, lsst::jointcal::SimpleAstrometryModel

class AstrometryTransform
#include <AstrometryTransform.h>

a virtual (interface) class for geometric transformations.

We implement here One AstrometryTransform interface class, and actual derived classes. Composition in the usual (mathematical) sense is provided using compose(), and some classes (e.g. AstrometryTransformLinear

) handle a * operator. Generic inversion by iteration exists, but it is at least 10 times slower than the corresponding “direct

transformation”. If a transform has an analytical inverse, then providing inverseTransform is obviously a very good idea. Before resorting to inverseTransform, consider using

StarMatchList::inverseTransform(). AstrometryTransformLinear::inverted() and TanPixelToRaDec::inverted() exist. The classes also provide derivation and linear approximation.

Subclassed by lsst::jointcal::AstrometryTransformIdentity, lsst::jointcal::AstrometryTransformPolynomial, lsst::jointcal::AstrometryTransformSkyWcs, lsst::jointcal::BaseTanWcs, lsst::jointcal::TanRaDecToPixel, lsst::jointcal::UserTransform

class AstrometryTransformIdentity : public lsst::jointcal::AstrometryTransform
#include <AstrometryTransform.h>

A do-nothing transformation. It anyway has dummy routines to mimick a AstrometryTransform.

class AstrometryTransformLinear : public lsst::jointcal::AstrometryTransformPolynomial
#include <AstrometryTransform.h>

implements the linear transformations (6 real coefficients).

Subclassed by lsst::jointcal::AstrometryTransformLinearRot, lsst::jointcal::AstrometryTransformLinearScale, lsst::jointcal::AstrometryTransformLinearShift

class AstrometryTransformLinearRot : public lsst::jointcal::AstrometryTransformLinear
#include <AstrometryTransform.h>

just here to provide a specialized constructor, and fit.

class AstrometryTransformLinearScale : public lsst::jointcal::AstrometryTransformLinear
#include <AstrometryTransform.h>

just here to provide specialized constructors. AstrometryTransformLinear fit routine.

class AstrometryTransformLinearShift : public lsst::jointcal::AstrometryTransformLinear
#include <AstrometryTransform.h>

just here to provide a specialized constructor, and fit.

class AstrometryTransformPolynomial : public lsst::jointcal::AstrometryTransform
#include <AstrometryTransform.h>

Polynomial transformation class.

Subclassed by lsst::jointcal::AstrometryTransformLinear

class AstrometryTransformSkyWcs : public lsst::jointcal::AstrometryTransform
#include <AstrometryTransform.h>

A AstrometryTransform that holds a SkyWcs

This is intended to hold the initial estimate for the WCS. It need not be TAN-SIP, nor exactly representable as a FITS WCS.

AstrometryTransformSkyWcs does not inherit from BaseTanWcs for two reasons:

  • There is no need.

  • It is not clear how to implement all of the BaseTanWcs interface, especially corrections.

class BaseStar : public lsst::jointcal::FatPoint
#include <BaseStar.h>

The base class for handling stars. Used by all matching routines.

Subclassed by lsst::jointcal::FittedStar, lsst::jointcal::MeasuredStar, lsst::jointcal::RefStar

class BaseTanWcs : public lsst::jointcal::AstrometryTransform

Subclassed by lsst::jointcal::TanPixelToRaDec, lsst::jointcal::TanSipPixelToRaDec

class CcdImage
#include <CcdImage.h>

Handler of an actual image from a single CCD. NOTE: could possibly be replaced with a subclass of afw.image.Exposure?

struct CcdImageKey
#include <CcdImage.h>

For hashing a ccdImage: the pair of (visit, ccd) IDs should be unique to each ccdImage.

class Chi2Accumulator
#include <Chi2.h>

Base class for Chi2Statistic and Chi2List, to allow addEntry inside Fitter for either class.

Essentially a mixin.

Subclassed by lsst::jointcal::Chi2List, lsst::jointcal::Chi2Statistic

class Chi2List : public lsst::jointcal::Chi2Accumulator, public std::vector<Chi2Star>
#include <Chi2.h>

Structure to accumulate the chi2 contributions per each star (to help find outliers).

class Chi2Statistic : public lsst::jointcal::Chi2Accumulator
#include <Chi2.h>

Simple structure to accumulate chi2 and ndof.

class ChipVisitAstrometryMapping : public lsst::jointcal::AstrometryMapping
#include <ChipVisitAstrometryMapping.h>

The mapping with two transforms in a row.

class ChipVisitPhotometryMapping : public lsst::jointcal::PhotometryMappingBase
#include <PhotometryMapping.h>

A two-level photometric transform: one for the ccd and one for the visit.

Subclassed by lsst::jointcal::ChipVisitFluxMapping, lsst::jointcal::ChipVisitMagnitudeMapping

class ConstrainedAstrometryModel : public lsst::jointcal::AstrometryModel
#include <ConstrainedAstrometryModel.h>

A multi-component model, fitting mappings for sensors and visits simultaneously.

This is the model used to fit mappings as the combination of a transformation depending on the chip number (instrument model) and a transformation per visit (anamorphism). The two-transformation Mapping required for this model is ChipVisitAstrometryMapping. This modeling of distortions is meant for a set of images from a single mosaic imager.

Parameters
  • ccdImageList: The exposures that will be fit.

  • projectionHandler: The projection from “Sky” (where the “true” coordinates live) to “Tangent Plane” (where the fitting occurs).

  • chipOrder: The polynomial order of the pixel->focal plane mapping for each sensor.

  • visitOrder: The polynomial order of the focal plane->tangent plane mapping for each visit.

class ConstrainedPhotometryModel : public lsst::jointcal::PhotometryModel
#include <ConstrainedPhotometryModel.h>

Photometry model with constraints, \(M(x,y) = M_CCD(x,y)*M_visit(u,v)\)

This model consists of the following components:

  • A spatially invariant zero point per CCD, constrained across all visits, \(M_CCD\).

  • A Chebyshev polynomial ( \(a_ij*T_i(x)*T_j(y)\) ) per visit, constrained across all CCDs, \(M_visit\).

Because this model’s parameters are degenerate under multiplication by a constant, \(M=(a*M_CCD)*(1/a*M_visit)\), we hold one CCD’s zero point fixed to remove that degeneracy.

Subclassed by lsst::jointcal::ConstrainedFluxModel, lsst::jointcal::ConstrainedMagnitudeModel

class FastFinder
#include <FastFinder.h>

Fast locator in starlists.

This is an auxillary class for matching objects from starlists. It allows to locate rapidly the closest objects from a given position. The very simple strategy is to sort objects according to 1 coordinate x, and to build an index that allows to select the objects with the x coordinate inside an interval. Then every slice in x is sorted according to y, which enables a fast scan inside a x slice. listMatchCollect takes about 10ms (PC 450 MHz, optimized “-O4”) for a match between lists of about 2000 objects each, which is fast enough for our needs. The same “locator” is used in listMatchupShift, to avoid scanning the whole input lists. Timing on listMatchCollect and listMatchupShift indicates a gain in speed by more than one order of magnitude after implementation of this FastFinder.

class FatPoint : public lsst::jointcal::Point
#include <FatPoint.h>

A Point with uncertainties.

Subclassed by lsst::jointcal::BaseStar

class FittedStar : public lsst::jointcal::BaseStar, public lsst::jointcal::PmBlock
#include <FittedStar.h>

The objects which have been measured several times.

MeasuredStars from different CcdImages that represent the same on-sky object all point to one FittedStar.

class FittedStarList : public lsst::jointcal::StarList<FittedStar>
#include <FittedStar.h>

A list of FittedStar s. Such a list is typically constructed by Associations.

class FitterBase
#include <FitterBase.h>

Base class for fitters.

Implements minimize and findOutliers. Chi2, residual, derivative, etc. calculations must be implemented in the child class via the virtual methods.

Subclassed by lsst::jointcal::AstrometryFit, lsst::jointcal::PhotometryFit

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

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

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

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

initialCalibFlux * SpatiallyInvariantTransform -> correctedFlux

class Frame
#include <Frame.h>

rectangle with sides parallel to axes.

when Frame’s are used to define subparts of images, xMin and xMax refer to the first and last pixels in the subimage

class IdentityProjectionHandler : public lsst::jointcal::ProjectionHandler
#include <ProjectionHandler.h>

The simplest implementation of ProjectionHandler. Means that coordinates of objects are expressed in the same space as the arrival mapping space. This is useful for fitting transform rms between images.

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

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

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

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

initialMagnitude + SpatiallyInvariantTransform -> correctedMagnitude

struct MatchConditions
#include <ListMatch.h>

Parameters to be provided to combinatorial searches.

class MeasuredStar : public lsst::jointcal::BaseStar
#include <MeasuredStar.h>

objects measured on actual images. Coordinates and uncertainties are expressed in pixel image frame. Flux expressed in ADU/s.

class MeasuredStarList : public lsst::jointcal::StarList<MeasuredStar>
#include <MeasuredStar.h>

A list of MeasuredStar. They are usually filled in Associations::createCcdImage.

class OneTPPerVisitHandler : public lsst::jointcal::ProjectionHandler
#include <ProjectionHandler.h>

A projection handler in which all CCDs from the same visit have the same tangent point.

We arbitrarily chose that all chips from a same visit have the same tangent point.

class PhotometryFit : public lsst::jointcal::FitterBase
#include <PhotometryFit.h>

Class that handles the photometric least squares problem.

class PhotometryMapping : public lsst::jointcal::PhotometryMappingBase
#include <PhotometryMapping.h>

A mapping containing a single photometryTransform.

class PhotometryMappingBase
#include <PhotometryMapping.h>

Relates transform(s) to their position in the fitting matrix and allows interaction with the transform(s).

Subclassed by lsst::jointcal::ChipVisitPhotometryMapping, lsst::jointcal::PhotometryMapping

class PhotometryModel

Subclassed by lsst::jointcal::ConstrainedPhotometryModel, lsst::jointcal::SimplePhotometryModel

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

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

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

Photometry offset independent of position. Abstract class.

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

struct PmBlock
#include <FittedStar.h>

objects whose position is going to be fitted. Coordinates in Common Tangent Plane.

Subclassed by lsst::jointcal::FittedStar

class Point
#include <Point.h>

A point in a plane.

Subclassed by lsst::jointcal::FatPoint

struct ProjectionHandler
#include <ProjectionHandler.h>

This is a virtual class that allows a lot of freedom in the choice of the projection from “Sky” (where coodinates are reported) to tangent plane (where they are compared to transformed measurements)

Subclassed by lsst::jointcal::IdentityProjectionHandler, lsst::jointcal::OneTPPerVisitHandler

class RefStar : public lsst::jointcal::BaseStar
#include <RefStar.h>

Objects used as position anchors, typically USNO stars. Coordinate system defined by user. The Common Tangent Plane seems a good idea.

class SimpleAstrometryMapping : public lsst::jointcal::AstrometryMapping

Subclassed by lsst::jointcal::SimplePolyMapping

class SimpleAstrometryModel : public lsst::jointcal::AstrometryModel
#include <SimpleAstrometryModel.h>

A model where there is one independent transform per CcdImage.

This modeling of distortions can even accommodate images set mixing instruments

Parameters
  • ccdImageList: The exposures that will be fit.

  • projectionHandler: The projection from “Sky” (where the “true” coordinates live) to “Tangent Plane” (where the fitting occurs).

  • initFromWCS: Initialize the model parameters from the original exposure Wcs parameters?

  • nNotFit: How many exposure to hold fixed and not be fit? (the first n will be selected) .

  • order: The polynomial order of each exposure’s pixel-tangent plane mapping.

class SimplePhotometryModel : public lsst::jointcal::PhotometryModel
#include <SimplePhotometryModel.h>

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

Subclassed by lsst::jointcal::SimpleFluxModel, lsst::jointcal::SimpleMagnitudeModel

class SimplePolyMapping : public lsst::jointcal::SimpleAstrometryMapping
#include <SimpleAstrometryMapping.h>

Mapping implementation for a polynomial transformation.

class SparseHisto4d
#include <Histo4d.h>

A class to histogram in 4 dimensions. Uses Sparse storage. The number of bin is limited to 256 per dimension. Used in ListMatch.cc

template<class Star>
class StarList : public std::list<std::shared_ptr<Star>>
#include <StarList.h>

std::lists of Stars.

It is a template class, which means that the Star type remains undefined until a user defines it. The std::list related operations (insertion, sort, traversal) are to be carried out using STL list operations. Most of the Star operations rely on routines to be provided in the Star class, usually user defined. The instanciation of this class for BaseStar (i.e. the replacement of the formal parameter ‘Star’ by ‘BaseStar’) is called BaseStarList. Take care: what is stored is pointers on Star’s and NOT Star’s. This implies that Stars being inserted in the std::list have to be obtained using ‘new’.

class StarMatch
#include <StarMatch.h>

A hanger for star associations.

class TanPixelToRaDec : public lsst::jointcal::BaseTanWcs
#include <AstrometryTransform.h>

The transformation that handles pixels to sideral transformations (Gnomonic, possibly with polynomial distortions).

class TanRaDecToPixel : public lsst::jointcal::AstrometryTransform
#include <AstrometryTransform.h>

This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane)

this transform does not implement corrections, since they are defined the other way around (from pixels to sky), and not invertible analytically. The inversion of tangent point WCS (TanPixelToRaDec) is obtained via inverseTransform().

class TanSipPixelToRaDec : public lsst::jointcal::BaseTanWcs
#include <AstrometryTransform.h>

Implements the (forward) SIP distorsion scheme.

class UserTransform : public lsst::jointcal::AstrometryTransform
#include <AstrometryTransform.h>

A run-time transform that allows users to define a AstrometryTransform with minimal coding (just the apply routine).